PLASMA
2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
|
This is the group of double complex functions using the advanced synchronous user interface.
int PLASMA_zcgels_Tile | ( | PLASMA_enum | trans, |
PLASMA_desc * | A, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | B, | ||
PLASMA_desc * | X, | ||
int * | ITER | ||
) |
PLASMA_zcgels_Tile - Solves overdetermined or underdetermined linear system of equations using the tile QR or the tile LQ factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_zcgesv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | trans | Intended usage: = PlasmaNoTrans: the linear system involves A; = PlasmaConjTrans: the linear system involves A**H. Currently only PlasmaNoTrans is supported. |
[in,out] | A |
|
[out] | T | On exit:
|
[in,out] | B | On entry, the M-by-NRHS matrix B of right hand side vectors, stored columnwise; On exit, if return value = 0, B is overwritten by the solution vectors, stored columnwise: if M >= N, rows 1 to N of B contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of the modulus of elements N+1 to M in that column; if M < N, rows 1 to N of B contain the minimum norm solution vectors; |
[in] | B | The descriptor of the M-by-NRHS matrix B of right hand side vectors, stored columnwise. Not modified. |
[in,out] | X | On entry, it's only the descriptor where to store the result. On exit, if return value = 0, X is the solution vectors, stored columnwise: if M >= N, rows 1 to N of X contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of the modulus of elements N+1 to M in that column; if M < N, rows 1 to N of X contain the minimum norm solution vectors; |
[out] | ITER | If > 0, ITER is the number of the current iteration in the iterative refinement process. -ITERMAX-1, if the refinment step failed and the double precision factorization has been used. |
PLASMA_SUCCESS | successful exit |
Definition at line 390 of file zcgels.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, PLASMA_zcgels_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zcgesv_Tile | ( | PLASMA_desc * | A, |
int * | IPIV, | ||
PLASMA_desc * | B, | ||
PLASMA_desc * | X, | ||
int * | ITER | ||
) |
PLASMA_zcgesv_Tile - Solves a system of linear equations using the tile LU factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_zcgesv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the N-by-N coefficient matrix A.
|
[out] | L | On exit:
|
[out] | IPIV | On exit, the pivot indices that define the permutations (not equivalent to LAPACK). |
[in] | B | On entry, the N-by-NRHS matrix of right hand side matrix B. |
[out] | X | On exit, if return value = 0, the N-by-NRHS solution matrix X. |
[out] | ITER | The number of the current iteration in the iterative refinement process |
PLASMA_SUCCESS | successful exit |
>0 | if i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed. |
Definition at line 374 of file zcgesv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, PLASMA_zcgesv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zcposv_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
PLASMA_desc * | X, | ||
int * | ITER | ||
) |
PLASMA_zcposv_Tile - Solves a symmetric positive definite or Hermitian positive definite system of linear equations using the Cholesky factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_zcposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in,out] | A | On entry, the N-by-N symmetric positive definite (or Hermitian) 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 = '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.
|
[in] | B | On entry, the N-by-NRHS matrix of right hand side matrix B. |
[out] | X | On exit, if return value = 0, the N-by-NRHS solution matrix X. |
[out] | ITER | The number of the current iteration in the iterative refinement process |
PLASMA_SUCCESS | successful exit |
>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. |
Definition at line 323 of file zcposv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, PLASMA_zcposv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zcungesv_Tile | ( | PLASMA_enum | trans, |
PLASMA_desc * | A, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | B, | ||
PLASMA_desc * | X, | ||
int * | ITER | ||
) |
PLASMA_zcungesv_Tile - Solves symmetric linear system of equations using the tile QR or the tile LQ factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_zcungesv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | trans | Intended usage: = PlasmaNoTrans: the linear system involves A; = PlasmaConjTrans: the linear system involves A**H. Currently only PlasmaNoTrans is supported. |
[in,out] | A |
|
[out] | T | On exit:
|
[in,out] | B | On entry, the M-by-NRHS matrix B of right hand side vectors, stored columnwise; |
[in] | B | The N-by-NRHS matrix of right hand side matrix B. |
[out] | X | If return value = 0, X is the solution vectors, stored columnwise: if M >= N, rows 1 to N of X contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of the modulus of elements N+1 to M in that column; if M < N, rows 1 to N of X contain the minimum norm solution vectors; |
[out] | ITER | The number of the current iteration in the iterative refinement process |
PLASMA_SUCCESS | successful exit |
Definition at line 328 of file zcungesv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, PLASMA_zcungesv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgebrd_Tile | ( | PLASMA_enum | jobu, |
PLASMA_enum | jobvt, | ||
PLASMA_desc * | A, | ||
double * | D, | ||
double * | E, | ||
PLASMA_desc * | U, | ||
PLASMA_desc * | VT, | ||
PLASMA_desc * | T | ||
) |
PLASMA_zgebrd_Tile - computes the singular value decomposition (SVD) of a complex M-by-N matrix A, optionally computing the left and/or right singular vectors. Tile equivalent of PLASMA_zgebrd(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | jobu | Specifies options for computing all or part of the matrix U. Intended usage: = PlasmaVec: all M columns of U are returned in array U; = PlasmaNoVec: no columns of U (no left singular vectors) are computed. |
[in] | jobvt | Specifies options for computing all or part of the matrix V**H. Intended usage: = PlasmaVec: all M columns of U are returned in array U; = PlasmaNoVec: no columns of U (no left singular vectors) are computed. |
[in,out] | A | On entry, the M-by-N matrix A. On exit, if JOBU = 'O', A is overwritten with the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBVT = 'O', A is overwritten with the first min(m,n) rows of V**H (the right singular vectors, stored rowwise); if JOBU .ne. 'O' and JOBVT .ne. 'O', the contents of A are destroyed. |
[out] | S | The singular values of A, sorted so that S(i) >= S(i+1). |
[out] | U | (LDU,M) if JOBU = 'A' or (LDU,min(M,N)) if JOBU = 'S'. If JOBU = 'A', U contains the M-by-M unitary matrix U; if JOBU = 'S', U contains the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBU = 'N' or 'O', U is not referenced. |
[out] | VT | If JOBVT = 'A', VT contains the N-by-N unitary matrix V**H; if JOBVT = 'S', VT contains the first min(m,n) rows of V**H (the right singular vectors, stored rowwise); if JOBVT = 'N' or 'O', VT is not referenced. |
[out] | T | On exit, auxiliary factorization data. |
Definition at line 333 of file zgebrd.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgebrd_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgelqf_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | T | ||
) |
PLASMA_zgelqf_Tile - Computes the tile LQ factorization of a matrix. Tile equivalent of PLASMA_zgelqf(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the M-by-N matrix A. On exit, the elements on and below the diagonal of the array contain the m-by-min(M,N) lower trapezoidal matrix L (L is lower triangular if M <= N); the elements above the diagonal represent the unitary matrix Q as a product of elementary reflectors, stored by tiles. |
[out] | T | On exit, auxiliary factorization data, required by PLASMA_zgelqs to solve the system of equations. |
PLASMA_SUCCESS | successful exit |
Definition at line 187 of file zgelqf.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgelqf_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgelqs_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | T, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgelqs_Tile - Computes a minimum-norm solution using previously computed LQ factorization. Tile equivalent of PLASMA_zgelqs(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | Details of the LQ factorization of the original matrix A as returned by PLASMA_zgelqf. |
[in] | T | Auxiliary factorization data, computed by PLASMA_zgelqf. |
[in,out] | B | On entry, the M-by-NRHS right hand side matrix B. On exit, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 207 of file zgelqs.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgelqs_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgels_Tile | ( | PLASMA_enum | trans, |
PLASMA_desc * | A, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgels_Tile - Solves overdetermined or underdetermined linear system of equations using the tile QR or the tile LQ factorization. Tile equivalent of PLASMA_zgels(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | trans | Intended usage: = PlasmaNoTrans: the linear system involves A; = PlasmaConjTrans: the linear system involves A**H. Currently only PlasmaNoTrans is supported. |
[in,out] | A | On entry, the M-by-N matrix A. On exit, if M >= N, A is overwritten by details of its QR factorization as returned by PLASMA_zgeqrf; if M < N, A is overwritten by details of its LQ factorization as returned by PLASMA_zgelqf. |
[out] | T | On exit, auxiliary factorization data. |
[in,out] | B | On entry, the M-by-NRHS matrix B of right hand side vectors, stored columnwise; On exit, if return value = 0, B is overwritten by the solution vectors, stored columnwise: if M >= N, rows 1 to N of B contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of the modulus of elements N+1 to M in that column; if M < N, rows 1 to N of B contain the minimum norm solution vectors; |
Definition at line 267 of file zgels.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgels_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgemm_Tile | ( | PLASMA_enum | transA, |
PLASMA_enum | transB, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
PLASMA_Complex64_t | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zgemm_Tile - Performs matrix multiplication. Tile equivalent of PLASMA_zgemm(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | transA | Specifies whether the matrix A is transposed, not transposed or conjugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is transposed; = PlasmaConjTrans: A is conjugate transposed. |
[in] | transB | Specifies whether the matrix B is transposed, not transposed or conjugate transposed: = PlasmaNoTrans: B is not transposed; = PlasmaTrans: B is transposed; = PlasmaConjTrans: B is conjugate transposed. |
[in] | alpha | alpha specifies the scalar alpha |
[in] | A | A is a LDA-by-ka matrix, where ka is K when transA = PlasmaNoTrans, and is M otherwise. |
[in] | B | B is a LDB-by-kb matrix, where kb is N when transB = PlasmaNoTrans, and is K otherwise. |
[in] | beta | beta specifies the scalar beta |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array is overwritten by the M by N matrix ( alpha*op( A )*op( B ) + beta*C ) |
PLASMA_SUCCESS | successful exit |
Definition at line 264 of file zgemm.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgemm_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgeqrf_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | T | ||
) |
PLASMA_zgeqrf_Tile - Computes the tile QR factorization of a matrix. Tile equivalent of PLASMA_zgeqrf(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N); the elements below the diagonal represent the unitary matrix Q as a product of elementary reflectors stored by tiles. |
[out] | T | On exit, auxiliary factorization data, required by PLASMA_zgeqrs to solve the system of equations. |
PLASMA_SUCCESS | successful exit |
Definition at line 186 of file zgeqrf.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgeqrf_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgeqrs_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | T, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgeqrs_Tile - Computes a minimum-norm solution using the tile QR factorization. Tile equivalent of PLASMA_zgeqrf(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | Details of the QR factorization of the original matrix A as returned by PLASMA_zgeqrf. |
[in] | T | Auxiliary factorization data, computed by PLASMA_zgeqrf. |
[in,out] | B | On entry, the m-by-nrhs right hand side matrix B. On exit, the n-by-nrhs solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 206 of file zgeqrs.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgeqrs_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgesv_incpiv_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | L, | ||
int * | IPIV, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgesv_incpiv_Tile - Solves a system of linear equations using the tile LU factorization. Tile equivalent of PLASMA_zgetrf_incpiv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the N-by-N coefficient matrix A. On exit, the tile L and U factors from the factorization (not equivalent to LAPACK). |
[in,out] | L | On exit, auxiliary factorization data, related to the tile L factor, necessary to solve the system of equations. |
[out] | IPIV | On exit, the pivot indices that define the permutations (not equivalent to LAPACK). |
[in,out] | B | On entry, the N-by-NRHS matrix of right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
>0 | if i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed. |
Definition at line 203 of file zgesv_incpiv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgesv_incpiv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgesv_Tile | ( | PLASMA_desc * | A, |
int * | IPIV, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgesv_Tile - Solves a system of linear equations using the tile LU factorization. Tile equivalent of PLASMA_zgetrf(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the N-by-N coefficient matrix A. On exit, the tile L and U factors from the factorization. |
[out] | IPIV | On exit, the pivot indices that define the permutations. |
[in,out] | B | On entry, the N-by-NRHS matrix of right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
>0 | if i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution could not be computed. |
Definition at line 187 of file zgesv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgesv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgesvd_Tile | ( | PLASMA_enum | jobu, |
PLASMA_enum | jobvt, | ||
PLASMA_desc * | A, | ||
double * | S, | ||
PLASMA_desc * | U, | ||
PLASMA_desc * | VT, | ||
PLASMA_desc * | T | ||
) |
PLASMA_zgesvd_Tile - computes the singular value decomposition (SVD) of a complex M-by-N matrix A, optionally computing the left and/or right singular vectors. Tile equivalent of PLASMA_zgesvd(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | jobu | Specifies options for computing all or part of the matrix U. Intended usage: = PlasmaVec: all M columns of U are returned in array U; = PlasmaNoVec: no columns of U (no left singular vectors) are computed. |
[in] | jobvt | Specifies options for computing all or part of the matrix V**H. Intended usage: = PlasmaVec: all M columns of U are returned in array U; = PlasmaNoVec: no columns of U (no left singular vectors) are computed. |
[in,out] | A | On entry, the M-by-N matrix A. On exit, if JOBU = 'O', A is overwritten with the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBVT = 'O', A is overwritten with the first min(m,n) rows of V**H (the right singular vectors, stored rowwise); if JOBU .ne. 'O' and JOBVT .ne. 'O', the contents of A are destroyed. |
[out] | S | The singular values of A, sorted so that S(i) >= S(i+1). |
[out] | U | (LDU,M) if JOBU = 'A' or (LDU,min(M,N)) if JOBU = 'S'. If JOBU = 'A', U contains the M-by-M unitary matrix U; if JOBU = 'S', U contains the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBU = 'N' or 'O', U is not referenced. |
[out] | VT | If JOBVT = 'A', VT contains the N-by-N unitary matrix V**H; if JOBVT = 'S', VT contains the first min(m,n) rows of V**H (the right singular vectors, stored rowwise); if JOBVT = 'N' or 'O', VT is not referenced. |
[out] | T | On exit, auxiliary factorization data. |
Definition at line 333 of file zgesvd.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgesvd_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgetrf_incpiv_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | L, | ||
int * | IPIV | ||
) |
PLASMA_zgetrf_incpiv_Tile - Computes the tile LU factorization of a matrix. Tile equivalent of PLASMA_zgetrf_incpiv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the M-by-N matrix to be factored. On exit, the tile factors L and U from the factorization. |
[out] | L | On exit, auxiliary factorization data, related to the tile L factor, required by PLASMA_zgetrs_incpiv to solve the system of equations. |
[out] | IPIV | The pivot indices that define the permutations (not equivalent to LAPACK). |
PLASMA_SUCCESS | successful exit |
>0 | if i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations. |
Definition at line 184 of file zgetrf_incpiv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgetrf_incpiv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgetrf_Tile | ( | PLASMA_desc * | A, |
int * | IPIV | ||
) |
PLASMA_zgetrf_Tile - Computes the tile LU factorization of a matrix. Tile equivalent of PLASMA_zgetrf(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the M-by-N matrix to be factored. On exit, the tile factors L and U from the factorization. |
[out] | L | On exit, auxiliary factorization data, related to the tile L factor, required by PLASMA_zgetrs to solve the system of equations. |
[out] | IPIV | The pivot indices that define the permutations (not equivalent to LAPACK). |
PLASMA_SUCCESS | successful exit |
>0 | if i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations. |
Definition at line 189 of file zgetrf.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgetrf_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgetri_Tile | ( | PLASMA_desc * | A, |
int * | IPIV | ||
) |
PLASMA_zgetri_Tile - Computes the inverse of a matrix using the LU factorization computed by PLASMA_zgetrf. This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A). Tile equivalent of PLASMA_zgetri(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in,out] | A | On entry, the triangular factor L or U from the factorization A = P*L*U as computed by PLASMA_zgetrf. On exit, if return value = 0, the inverse of the original matrix A. |
[in] | IPIV | The pivot indices that define the permutations as returned by PLASMA_zgetrf. |
PLASMA_SUCCESS | successful exit |
>0 | if i, the (i,i) element of the factor U is exactly zero; The matrix is singular and its inverse could not be computed. |
Definition at line 175 of file zgetri.c.
References PLASMA_Alloc_Workspace_zgetri_Tile_Async(), plasma_context_self(), plasma_desc_mat_free(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgetri_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgetrs_incpiv_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | L, | ||
int * | IPIV, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgetrs_incpiv_Tile - Solves a system of linear equations using previously computed LU factorization. Tile equivalent of PLASMA_zgetrs_incpiv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | The tile factors L and U from the factorization, computed by PLASMA_zgetrf_incpiv. |
[in] | L | Auxiliary factorization data, related to the tile L factor, computed by PLASMA_zgetrf_incpiv. |
[in] | IPIV | The pivot indices from PLASMA_zgetrf_incpiv (not equivalent to LAPACK). |
[in,out] | B | On entry, the N-by-NRHS matrix of right hand side matrix B. On exit, the solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 206 of file zgetrs_incpiv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgetrs_incpiv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zgetrs_Tile | ( | PLASMA_enum | trans, |
PLASMA_desc * | A, | ||
int * | IPIV, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zgetrs_Tile - Solves a system of linear equations using previously computed LU factorization. Tile equivalent of PLASMA_zgetrs(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | trans | Intended to specify the the form of the system of equations: = PlasmaNoTrans: A * X = B (No transpose) = PlasmaTrans: A**T * X = B (Transpose) = PlasmaConjTrans: A**H * X = B (Conjugate transpose) |
[in] | A | The tile factors L and U from the factorization, computed by PLASMA_zgetrf. |
[in] | IPIV | The pivot indices from PLASMA_zgetrf. |
[in,out] | B | On entry, the N-by-NRHS matrix of right hand side matrix B. On exit, the solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 199 of file zgetrs.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zgetrs_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zheev_Tile | ( | PLASMA_enum | jobz, |
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
double * | W, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | Q | ||
) |
PLASMA_zheev_Tile - Computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A using a two-stage approach: First stage: reduction to band tridiagonal form; Second stage: reduction from band to tridiagonal form.
Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | jobz | Intended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. PlasmaVec is NOT supported. |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in,out] | A | On entry, the symmetric (or Hermitian) 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 jobz = PlasmaVec, then if return value = 0, A contains the orthonormal eigenvectors of the matrix A. If jobz = PlasmaNoVec, then on exit the lower triangle (if uplo = PlasmaLower) or the upper triangle (if uplo = PlasmaUpper) of A, including the diagonal, is destroyed.* |
[out] | W | On exit, if info = 0, the eigenvalues. |
[in,out] | T | On entry, descriptor as return by PLASMA_Alloc_Workspace_zheev On exit, contains auxiliary factorization data. |
[out] | Q | On exit, if jobz = PlasmaVec and info = 0, the eigenvectors. |
PLASMA_SUCCESS | successful exit |
<0 | if -i, the i-th argument had an illegal value |
>0 | if INFO = i, the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero. |
Definition at line 274 of file zheev.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zheev_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zheevd_Tile | ( | PLASMA_enum | jobz, |
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
double * | W, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | Q | ||
) |
PLASMA_zheevd_Tile - Computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A using a two-stage approach: First stage: reduction to band tridiagonal form; Second stage: reduction from band to tridiagonal form.
Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | jobz | Intended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. PlasmaVec is NOT supported. |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in,out] | A | On entry, the symmetric (or Hermitian) 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 jobz = PlasmaVec, then if return value = 0, A contains the orthonormal eigenvectors of the matrix A. If jobz = PlasmaNoVec, then on exit the lower triangle (if uplo = PlasmaLower) or the upper triangle (if uplo = PlasmaUpper) of A, including the diagonal, is destroyed.* |
[out] | T | On exit, auxiliary factorization data. |
[out] | W | On exit, if info = 0, the eigenvalues. |
[out] | Q | On exit, if jobz = PlasmaVec and info = 0, the eigenvectors. |
PLASMA_SUCCESS | successful exit |
<0 | if -i, the i-th argument had an illegal value |
>0 | if INFO = i, the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero. |
Definition at line 277 of file zheevd.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zheevd_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zhegst_Tile | ( | PLASMA_enum | itype, |
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zhegst_Tile - reduces a complex Hermitian-definite generalized eigenproblem to standard form. If PlasmaItype == 1, the problem is A*x = lambda*B*x, and A is overwritten by inv(U**H)*A*inv(U) or inv(L)*A*inv(L**H) If PlasmaItype == 2 or 3, the problem is A*B*x = lambda*x or B*A*x = lambda*x, and A is overwritten by U*A*U**H or L**H*A*L. B must have been previously factorized as U**H*U or L*L**H by PLASMA_ZPOTRF. ONLY PlasmaItype == 1 and PlasmaLower supported! Tile equivalent of PLASMA_zhegst(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | PlasmaItype | Intended usage: = 1: A*x=(lambda)*B*x = 2: A*Bx=(lambda)*x = 3: B*A*x=(lambda)*x Currently only PlasmaItype == 1 is supported. |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. Currently only PlasmaLower is supported. |
[in,out] | A | On entry, the symmetric (or Hermitian) 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, if return value == 0, the transformed matrix, stored in the same format as A. |
[in,out] | B | On entry, the triangular factor from the Cholesky factorization of B, as returned by PLASMA_ZPOTRF. |
PLASMA_SUCCESS | successful exit |
<0 | if -i, the i-th argument had an illegal value |
Definition at line 233 of file zhegst.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zhegst_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zhegv_Tile | ( | PLASMA_enum | itype, |
PLASMA_enum | jobz, | ||
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
double * | W, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | Q | ||
) |
PLASMA_zhegv_Tile - Computes all eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian-definite eigenproblem of the form: A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. Here A and B are assumed to be Hermitian and B is also positive definite.
Tile equivalent of PLASMA_zhegv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | PlasmaItype | Intended usage: = 1: A*x=(lambda)*B*x = 2: A*Bx=(lambda)*x = 3: B*A*x=(lambda)*x |
[in] | jobz | Intended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. Currently only PlasmaVec is NOT supported. |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A and B are stored; = PlasmaLower: Lower triangle of A and B are stored. |
[in] | N | The order of the matrix A. N >= 0. |
[in,out] | A | On entry, the symmetric (or Hermitian) 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, if jobz = PlasmaVec, then if return value = 0, A contains the matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**H*B*Z = I; if ITYPE = 3, Z**H*inv(B)*Z = I. If jobz = PlasmaNoVec, then on exit the lower triangle (if uplo = PlasmaLower) or the upper triangle (if uplo = PlasmaUpper) of A, including the diagonal, is destroyed. |
[in] | LDA | The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | On entry, the symmetric (or Hermitian) positive definite matrix B. If uplo = PlasmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B, and the strictly lower triangular part of B is not referenced. If uplo = PlasmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B, and the strictly upper triangular part of B is not referenced. On exit, if return value <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | LDB | The leading dimension of the array B. LDA >= max(1,N). |
[in,out] | T | On entry, descriptor as return by PLASMA_Alloc_Workspace_zhegv On exit, contains auxiliary factorization data. |
[out] | W | On exit, if info = 0, the eigenvalues. |
[out] | Q | On exit, if jobz = PlasmaVec and info = 0, the eigenvectors. |
PLASMA_SUCCESS | successful exit |
<0 | if -i, the i-th argument had an illegal value |
<=N | if INFO = i, plasma_zhegv failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero. |
>N | if INFO = N + i, for 1 <= i <= N, then the leading minor of order i of B is not positive definite. The factorization of B could not be completed and no eigenvalues or eigenvectors were computed. |
Definition at line 369 of file zhegv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zhegv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zhemm_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | uplo, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
PLASMA_Complex64_t | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zhemm_Tile - Performs Hermitian matrix multiplication. Tile equivalent of PLASMA_zhemm(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Specifies whether the hermitian matrix A appears on the left or right in the operation as follows: = PlasmaLeft: = PlasmaRight:
|
[in] | uplo | Specifies whether the upper or lower triangular part of the hermitian matrix A is to be referenced as follows: = PlasmaLower: Only the lower triangular part of the hermitian matrix A is to be referenced. = PlasmaUpper: Only the upper triangular part of the hermitian matrix A is to be referenced. |
[in] | alpha | Specifies the scalar alpha. |
[in] | A | A is a LDA-by-ka matrix, where ka is M when side = PlasmaLeft, and is N otherwise. Only the uplo triangular part is referenced. |
[in] | B | B is a LDB-by-N matrix, where the leading M-by-N part of the array B must contain the matrix B. |
[in] | beta | Specifies the scalar beta. |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array is overwritten by the M by N updated matrix. |
PLASMA_SUCCESS | successful exit |
Definition at line 254 of file zhemm.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zhemm_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zher2k_Tile | ( | PLASMA_enum | uplo, |
PLASMA_enum | trans, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
double | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zher2k_Tile - Performs hermitian rank k update. Tile equivalent of PLASMA_zher2k(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored. |
[in] | trans | Specifies whether the matrix A is transposed or conjugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaConjTrans: A is conjugate transposed. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | A is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise. |
[in] | B | B is a LDB-by-kb matrix, where kb is K when trans = PlasmaNoTrans, and is N otherwise. |
[in] | beta | beta specifies the scalar beta |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array uplo part of the matrix is overwritten by the uplo part of the updated matrix. |
PLASMA_SUCCESS | successful exit |
Definition at line 250 of file zher2k.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zher2k_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zherk_Tile | ( | PLASMA_enum | uplo, |
PLASMA_enum | trans, | ||
double | alpha, | ||
PLASMA_desc * | A, | ||
double | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zherk_Tile - Performs hermitian rank k update. Tile equivalent of PLASMA_zherk(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored. |
[in] | trans | Specifies whether the matrix A is transposed or conjugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaConjTrans: A is conjugate transposed. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | A is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise. |
[in] | beta | beta specifies the scalar beta |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array uplo part of the matrix is overwritten by the uplo part of the updated matrix. |
PLASMA_SUCCESS | successful exit |
Definition at line 227 of file zherk.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zherk_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zhetrd_Tile | ( | PLASMA_enum | jobz, |
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
double * | D, | ||
double * | E, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | Q | ||
) |
PLASMA_zhetrd_Tile - reduces a complex Hermitian matrix A to real symmetric tridiagonal form S using a two-stage approach First stage: reduction to band tridiagonal form (unitary Q1); Second stage: reduction from band to tridiagonal form (unitary Q2). Let Q = Q1 * Q2 be the global unitary transformation; Q**H * A * Q = S. Tile equivalent of PLASMA_zhetrd(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | jobz | Intended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. PlasmaVec is NOT supported. |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in,out] | A | On entry, the symmetric (or Hermitian) 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 jobz = PlasmaVec, then if return value = 0, A contains the orthonormal eigenvectors of the matrix A. If jobz = PlasmaNoVec, then on exit the lower triangle (if uplo = PlasmaLower) or the upper triangle (if uplo = PlasmaUpper) of A, including the diagonal, is destroyed.* |
[out] | D | On exit, the diagonal elements of the tridiagonal matrix: D(i) = A(i,i). |
[out] | E | On exit, he off-diagonal elements of the tridiagonal matrix: E(i) = A(i,i+1) if uplo = PlasmaUpper, E(i) = A(i+1,i) if uplo = PlasmaLower. |
[out] | T | On exit, auxiliary factorization data. |
[out] | Q | On exit, if jobz = PlasmaVec and info = 0, the eigenvectors. |
PLASMA_SUCCESS | successful exit |
<0 | if -i, the i-th argument had an illegal value |
>0 | if INFO = i, the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero. |
Definition at line 279 of file zhetrd.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zhetrd_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zlacpy_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zlacpy_Tile - Tile equivalent of PLASMA_zlacpy(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | Specifies the part of the matrix A to be copied to B. = PlasmaUpperLower: All the matrix A = PlasmaUpper: Upper triangular part = PlasmaLower: Lower triangular part |
[in] | A | The M-by-N matrix A. If uplo = PlasmaUpper, only the upper trapezium is accessed; if UPLO = PlasmaLower, only the lower trapezium is accessed. |
[out] | B | The M-by-N matrix B. On exit, B = A in the locations specified by UPLO. |
PLASMA_SUCCESS | successful exit |
Definition at line 183 of file zlacpy.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, and PLASMA_zlacpy_Tile_Async().
double PLASMA_zlange_Tile | ( | PLASMA_enum | norm, |
PLASMA_desc * | A, | ||
double * | work | ||
) |
PLASMA_zlange_Tile - Tile equivalent of PLASMA_zlange(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | norm | = PlasmaMaxNorm: Max norm = PlasmaOneNorm: One norm = PlasmaInfNorm: Infinity norm = PlasmaFrobeniusNorm: Frobenius norm |
[in] | A | On entry, the triangular factor U or L. On exit, if UPLO = 'U', the upper triangle of A is overwritten with the upper triangle of the product U * U'; if UPLO = 'L', the lower triangle of A is overwritten with the lower triangle of the product L' * L. |
[in] | work | double precision array of dimension (MAX(1,LWORK)), where LWORK >= M when NORM = PlasmaInfNorm; otherwise, WORK is not referenced. |
PLASMA_SUCCESS | successful exit |
Definition at line 193 of file zlange.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and PLASMA_zlange_Tile_Async().
double PLASMA_zlanhe_Tile | ( | PLASMA_enum | norm, |
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
double * | work | ||
) |
PLASMA_zlanhe_Tile - Tile equivalent of PLASMA_zlanhe(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | norm | = PlasmaMaxNorm: Max norm = PlasmaOneNorm: One norm = PlasmaInfNorm: Infinity norm = PlasmaFrobeniusNorm: Frobenius norm |
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | A | On entry, the triangular factor U or L. On exit, if UPLO = 'U', the upper triangle of A is overwritten with the upper triangle of the product U * U'; if UPLO = 'L', the lower triangle of A is overwritten with the lower triangle of the product L' * L. |
[in] | work | double precision array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise. |
PLASMA_SUCCESS | successful exit |
Definition at line 195 of file zlanhe.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and PLASMA_zlanhe_Tile_Async().
double PLASMA_zlansy_Tile | ( | PLASMA_enum | norm, |
PLASMA_enum | uplo, | ||
PLASMA_desc * | A, | ||
double * | work | ||
) |
PLASMA_zlansy_Tile - Tile equivalent of PLASMA_zlansy(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | norm | = PlasmaMaxNorm: Max norm = PlasmaOneNorm: One norm = PlasmaInfNorm: Infinity norm = PlasmaFrobeniusNorm: Frobenius norm |
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | A | On entry, the triangular factor U or L. On exit, if UPLO = 'U', the upper triangle of A is overwritten with the upper triangle of the product U * U'; if UPLO = 'L', the lower triangle of A is overwritten with the lower triangle of the product L' * L. |
[in] | work | double precision array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise. |
PLASMA_SUCCESS | successful exit |
Definition at line 195 of file zlansy.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and PLASMA_zlansy_Tile_Async().
int PLASMA_zlaset_Tile | ( | PLASMA_enum | uplo, |
PLASMA_Complex64_t | alpha, | ||
PLASMA_Complex64_t | beta, | ||
PLASMA_desc * | A | ||
) |
PLASMA_zlaset_Tile - Tile equivalent of PLASMA_zlaset(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | Specifies the part of the matrix A to be copied to B. = PlasmaUpperLower: All the matrix A = PlasmaUpper: Upper triangular part = PlasmaLower: Lower triangular part |
[in,out] | A | On entry, the m by n matrix A. On exit, A(i,j) = ALPHA, 1 <= i <= m, 1 <= j <= n, i.ne.j; A(i,i) = BETA , 1 <= i <= min(m,n) |
PLASMA_SUCCESS | successful exit |
Definition at line 172 of file zlaset.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, and PLASMA_zlaset_Tile_Async().
int PLASMA_zlaswp_Tile | ( | PLASMA_desc * | A, |
int | K1, | ||
int | K2, | ||
int * | IPIV, | ||
int | INCX | ||
) |
PLASMA_zlaswp_Tile - performs a series of row interchanges on the matrix A. One row interchange is initiated for each of rows K1 through K2 of A. Tile equivalent of PLASMA_zlaswp(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | The tile factors L and U from the factorization, computed by PLASMA_zgetrf. |
[in] | K1 | The first element of IPIV for which a row interchange will be done. |
[in] | K2 | The last element of IPIV for which a row interchange will be done. |
[in] | IPIV | The pivot indices from PLASMA_zgetrf. |
[in] | INCX | The increment between successive values of IPIV. If IPIV is negative, the pivots are applied in reverse order. |
PLASMA_SUCCESS | successful exit |
Definition at line 176 of file zlaswp.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zlaswp_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zlaswpc_Tile | ( | PLASMA_desc * | A, |
int | K1, | ||
int | K2, | ||
int * | IPIV, | ||
int | INCX | ||
) |
PLASMA_zlaswpc_Tile - performs a series of row interchanges on the matrix A. One row interchange is initiated for each of rows K1 through K2 of A. Tile equivalent of PLASMA_zlaswpc(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | The tile factors L and U from the factorization, computed by PLASMA_zgetrf. |
[in] | K1 | The first element of IPIV for which a row interchange will be done. |
[in] | K2 | The last element of IPIV for which a row interchange will be done. |
[in] | IPIV | The pivot indices from PLASMA_zgetrf. |
[in] | INCX | The increment between successive values of IPIV. If IPIV is negative, the pivots are applied in reverse order. |
PLASMA_SUCCESS | successful exit |
Definition at line 176 of file zlaswpc.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zlaswpc_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zlauum_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A | ||
) |
PLASMA_zlauum_Tile - Computes the product U * U' or L' * L, where the triangular factor U or L is stored in the upper or lower triangular part of the array A. Tile equivalent of PLASMA_zlauum(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | A | On entry, the triangular factor U or L. On exit, if UPLO = 'U', the upper triangle of A is overwritten with the upper triangle of the product U * U'; if UPLO = 'L', the lower triangle of A is overwritten with the lower triangle of the product L' * L. |
PLASMA_SUCCESS | successful exit |
Definition at line 172 of file zlauum.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zlauum_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zplghe_Tile | ( | double | bump, |
PLASMA_desc * | A, | ||
unsigned long long int | seed | ||
) |
PLASMA_zplghe_Tile - Generate a random hermitian matrix by tiles. Tile equivalent of PLASMA_zplghe(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | bump | The value to add to the diagonal to be sure to have a positive definite matrix. |
[in] | A | On exit, The random hermitian matrix A generated. |
[in] | seed | The seed used in the random generation. |
PLASMA_SUCCESS | successful exit |
Definition at line 153 of file zplghe.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zplghe_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zplgsy_Tile | ( | PLASMA_Complex64_t | bump, |
PLASMA_desc * | A, | ||
unsigned long long int | seed | ||
) |
PLASMA_zplgsy_Tile - Generate a random hermitian matrix by tiles. Tile equivalent of PLASMA_zplgsy(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | bump | The value to add to the diagonal to be sure to have a positive definite matrix. |
[in] | A | On exit, The random hermitian matrix A generated. |
[in] | seed | The seed used in the random generation. |
PLASMA_SUCCESS | successful exit |
Definition at line 153 of file zplgsy.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zplgsy_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zplrnt_Tile | ( | PLASMA_desc * | A, |
unsigned long long int | seed | ||
) |
PLASMA_zplrnt_Tile - Generate a random matrix by tiles. Tile equivalent of PLASMA_zplrnt(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | On exit, The random matrix A generated. |
[in] | seed | The seed used in the random generation. |
PLASMA_SUCCESS | successful exit |
Definition at line 151 of file zplrnt.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zplrnt_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zposv_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zposv_Tile - Solves a symmetric positive definite or Hermitian positive definite system of linear equations using the Cholesky factorization. Tile equivalent of PLASMA_zposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in,out] | A | On entry, the symmetric positive definite (or Hermitian) 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] | B | On entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
>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. |
Definition at line 212 of file zposv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zposv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zpotrf_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A | ||
) |
PLASMA_zpotrf_Tile - Computes the Cholesky factorization of a symmetric positive definite or Hermitian positive definite matrix. Tile equivalent of PLASMA_zpotrf(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | A | On entry, the symmetric positive definite (or Hermitian) 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. |
PLASMA_SUCCESS | successful exit |
>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. |
Definition at line 183 of file zpotrf.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zpotrf_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zpotri_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A | ||
) |
PLASMA_zpotri_Tile - 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 PLASMA_zpotrf. Tile equivalent of PLASMA_zpotri(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | A | On entry, the triangular factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H, as computed by PLASMA_zpotrf. On exit, the upper or lower triangle of the (Hermitian) inverse of A, overwriting the input factor U or L. |
PLASMA_SUCCESS | successful exit |
>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. |
Definition at line 172 of file zpotri.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zpotri_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zpotrs_Tile | ( | PLASMA_enum | uplo, |
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zpotrs_Tile - Solves a system of linear equations using previously computed Cholesky factorization. Tile equivalent of PLASMA_zpotrs(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | A | The triangular factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H, computed by PLASMA_zpotrf. |
[in,out] | B | On entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 187 of file zpotrs.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zpotrs_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zsymm_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | uplo, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
PLASMA_Complex64_t | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zsymm_Tile - Performs symmetric matrix multiplication. Tile equivalent of PLASMA_zsymm(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Specifies whether the symmetric matrix A appears on the left or right in the operation as follows: = PlasmaLeft: = PlasmaRight:
|
[in] | uplo | Specifies whether the upper or lower triangular part of the symmetric matrix A is to be referenced as follows: = PlasmaLower: Only the lower triangular part of the symmetric matrix A is to be referenced. = PlasmaUpper: Only the upper triangular part of the symmetric matrix A is to be referenced. |
[in] | alpha | Specifies the scalar alpha. |
[in] | A | A is a LDA-by-ka matrix, where ka is M when side = PlasmaLeft, and is N otherwise. Only the uplo triangular part is referenced. |
[in] | B | B is a LDB-by-N matrix, where the leading M-by-N part of the array B must contain the matrix B. |
[in] | beta | Specifies the scalar beta. |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array is overwritten by the M by N updated matrix. |
PLASMA_SUCCESS | successful exit |
Definition at line 254 of file zsymm.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zsymm_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zsyr2k_Tile | ( | PLASMA_enum | uplo, |
PLASMA_enum | trans, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B, | ||
PLASMA_Complex64_t | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zsyr2k_Tile - Performs symmetric rank k update. Tile equivalent of PLASMA_zsyr2k(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored. |
[in] | trans | Specifies whether the matrix A is transposed or conjugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is conjugate transposed. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | A is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise. |
[in] | B | B is a LDB-by-kb matrix, where kb is K when trans = PlasmaNoTrans, and is N otherwise. |
[in] | beta | beta specifies the scalar beta |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array uplo part of the matrix is overwritten by the uplo part of the updated matrix. |
PLASMA_SUCCESS | successful exit |
Definition at line 250 of file zsyr2k.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zsyr2k_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zsyrk_Tile | ( | PLASMA_enum | uplo, |
PLASMA_enum | trans, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_Complex64_t | beta, | ||
PLASMA_desc * | C | ||
) |
PLASMA_zsyrk_Tile - Performs rank k update. Tile equivalent of PLASMA_zsyrk(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored. |
[in] | trans | Specifies whether the matrix A is transposed or conjugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is transposed. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | A is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise. |
[in] | beta | beta specifies the scalar beta |
[in,out] | C | C is a LDC-by-N matrix. On exit, the array uplo part of the matrix is overwritten by the uplo part of the updated matrix. |
PLASMA_SUCCESS | successful exit |
Definition at line 227 of file zsyrk.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zsyrk_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_ztrmm_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | uplo, | ||
PLASMA_enum | transA, | ||
PLASMA_enum | diag, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_ztrmm_Tile - Computes triangular solve. Tile equivalent of PLASMA_ztrmm(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Specifies whether A appears on the left or on the right of X: = PlasmaLeft: A*X = B = PlasmaRight: X*A = B |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | transA | Specifies whether the matrix A is transposed, not transposed or conjugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaConjTrans: A is conjugate transposed. |
[in] | diag | Specifies whether or not A is unit triangular: = PlasmaNonUnit: A is non unit; = PlasmaUnit: A us unit. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | The triangular matrix A. If uplo = PlasmaUpper, the leading N-by-N upper triangular part of the array A contains the upper triangular matrix, and the strictly lower triangular part of A is not referenced. If uplo = PlasmaLower, the leading N-by-N lower triangular part of the array A contains the lower triangular matrix, and the strictly upper triangular part of A is not referenced. If diag = PlasmaUnit, the diagonal elements of A are also not referenced and are assumed to be 1. |
[in,out] | B | On entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 249 of file ztrmm.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_ztrmm_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_ztrsm_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | uplo, | ||
PLASMA_enum | transA, | ||
PLASMA_enum | diag, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_ztrsm_Tile - Computes triangular solve. Tile equivalent of PLASMA_ztrsm(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Specifies whether A appears on the left or on the right of X: = PlasmaLeft: A*X = B = PlasmaRight: X*A = B |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | transA | Specifies whether the matrix A is transposed, not transposed or conjugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaConjTrans: A is conjugate transposed. |
[in] | diag | Specifies whether or not A is unit triangular: = PlasmaNonUnit: A is non unit; = PlasmaUnit: A us unit. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | The triangular matrix A. If uplo = PlasmaUpper, the leading N-by-N upper triangular part of the array A contains the upper triangular matrix, and the strictly lower triangular part of A is not referenced. If uplo = PlasmaLower, the leading N-by-N lower triangular part of the array A contains the lower triangular matrix, and the strictly upper triangular part of A is not referenced. If diag = PlasmaUnit, the diagonal elements of A are also not referenced and are assumed to be 1. |
[in,out] | B | On entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 249 of file ztrsm.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_ztrsm_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_ztrsmpl_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | L, | ||
int * | IPIV, | ||
PLASMA_desc * | B | ||
) |
PLASMA_ztrsmpl_Tile - Performs the forward substitution step of solving a system of linear equations after the tile LU factorization of the matrix. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | The tile factor L from the factorization, computed by PLASMA_zgetrf_incpiv. |
[in] | L | Auxiliary factorization data, related to the tile L factor, computed by PLASMA_zgetrf_incpiv. |
[in] | IPIV | The pivot indices from PLASMA_zgetrf_incpiv (not equivalent to LAPACK). |
[in,out] | B | On entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 191 of file ztrsmpl.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_ztrsmpl_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_ztrsmrv_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | uplo, | ||
PLASMA_enum | transA, | ||
PLASMA_enum | diag, | ||
PLASMA_Complex64_t | alpha, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | B | ||
) |
PLASMA_ztrsmrv_Tile - Computes triangular solve. Tile equivalent of PLASMA_ztrsmrv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Specifies whether A appears on the left or on the right of X: = PlasmaLeft: A*X = B = PlasmaRight: X*A = B |
[in] | uplo | Specifies whether the matrix A is upper triangular or lower triangular: = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | transA | Specifies whether the matrix A is transposed, not transposed or conjugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaConjTrans: A is conjugate transposed. |
[in] | diag | Specifies whether or not A is unit triangular: = PlasmaNonUnit: A is non unit; = PlasmaUnit: A us unit. |
[in] | alpha | alpha specifies the scalar alpha. |
[in] | A | The triangular matrix A. If uplo = PlasmaUpper, the leading N-by-N upper triangular part of the array A contains the upper triangular matrix, and the strictly lower triangular part of A is not referenced. If uplo = PlasmaLower, the leading N-by-N lower triangular part of the array A contains the lower triangular matrix, and the strictly upper triangular part of A is not referenced. If diag = PlasmaUnit, the diagonal elements of A are also not referenced and are assumed to be 1. |
[in,out] | B | On entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X. |
PLASMA_SUCCESS | successful exit |
Definition at line 249 of file ztrsmrv.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_ztrsmrv_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_ztrtri_Tile | ( | PLASMA_enum | uplo, |
PLASMA_enum | diag, | ||
PLASMA_desc * | A | ||
) |
PLASMA_ztrtri_Tile - Computes the inverse of a complex upper or lower triangular matrix A. Tile equivalent of PLASMA_ztrtri(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | uplo | = PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored. |
[in] | diag | = PlasmaNonUnit: A is non-unit triangular; = PlasmaUnit: A us unit triangular. |
[in] | A | On entry, the triangular matrix A. If UPLO = 'U', the leading N-by-N upper triangular part of the array A contains the upper triangular matrix, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of the array A contains the lower triangular matrix, and the strictly upper triangular part of A is not referenced. If DIAG = 'U', the diagonal elements of A are also not referenced and are assumed to be 1. On exit, the (triangular) inverse of the original matrix, in the same storage format. |
PLASMA_SUCCESS | successful exit |
>0 | if i, A(i,i) is exactly zero. The triangular matrix is singular and its inverse can not be computed. |
Definition at line 191 of file ztrtri.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_ztrtri_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zunglq_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | T, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zunglq_Tile - Generates an M-by-N matrix Q with orthonormal rows, which is defined as the first M rows of a product of the elementary reflectors returned by PLASMA_zgelqf. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | Details of the LQ factorization of the original matrix A as returned by PLASMA_zgelqf. |
[in] | T | Auxiliary factorization data, computed by PLASMA_zgelqf. |
[out] | B | On exit, the M-by-N matrix Q. |
PLASMA_SUCCESS | successful exit |
Definition at line 202 of file zunglq.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zunglq_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zungqr_Tile | ( | PLASMA_desc * | A, |
PLASMA_desc * | T, | ||
PLASMA_desc * | Q | ||
) |
PLASMA_zungqr_Tile - Generates an M-by-N matrix Q with orthonormal columns, which is defined as the first N columns of a product of the elementary reflectors returned by PLASMA_zgeqrf. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | A | Details of the QR factorization of the original matrix A as returned by PLASMA_zgeqrf. |
[in] | T | Auxiliary factorization data, computed by PLASMA_zgeqrf. |
[out] | Q | On exit, the M-by-N matrix Q. |
PLASMA_SUCCESS | successful exit |
Definition at line 200 of file zungqr.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zungqr_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zunmlq_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | trans, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zunmlq_Tile - overwrites the general M-by-N matrix C with Q*C, where Q is an orthogonal matrix (unitary in the complex case) defined as the product of elementary reflectors returned by PLASMA_zgelqf_Tile Q is of order M. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Intended usage: = PlasmaLeft: apply Q or Q**H from the left; = PlasmaRight: apply Q or Q**H from the right. Currently only PlasmaLeft is supported. |
[in] | trans | Intended usage: = PlasmaNoTrans: no transpose, apply Q; = PlasmaConjTrans: conjugate transpose, apply Q**H. Currently only PlasmaConjTrans is supported. |
[in] | A | Details of the LQ factorization of the original matrix A as returned by PLASMA_zgelqf. |
[in] | T | Auxiliary factorization data, computed by PLASMA_zgelqf. |
[in,out] | B | On entry, the M-by-N matrix B. On exit, B is overwritten by Q*B or Q**H*B. |
PLASMA_SUCCESS | successful exit |
Definition at line 247 of file zunmlq.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zunmlq_Tile_Async(), and plasma_sequence_t::status.
int PLASMA_zunmqr_Tile | ( | PLASMA_enum | side, |
PLASMA_enum | trans, | ||
PLASMA_desc * | A, | ||
PLASMA_desc * | T, | ||
PLASMA_desc * | B | ||
) |
PLASMA_zunmqr_Tile - overwrites the general M-by-N matrix C with Q*C, where Q is an orthogonal matrix (unitary in the complex case) defined as the product of elementary reflectors returned by PLASMA_zgeqrf_Tile Q is of order M. All matrices are passed through descriptors. All dimensions are taken from the descriptors.
[in] | side | Intended usage: = PlasmaLeft: apply Q or Q**H from the left; = PlasmaRight: apply Q or Q**H from the right. Currently only PlasmaLeft is supported. |
[in] | trans | Intended usage: = PlasmaNoTrans: no transpose, apply Q; = PlasmaConjTrans: conjugate transpose, apply Q**H. Currently only PlasmaConjTrans is supported. |
[in] | A | Details of the QR factorization of the original matrix A as returned by PLASMA_zgeqrf. |
[in] | T | Auxiliary factorization data, computed by PLASMA_zgeqrf. |
[in,out] | B | On entry, the M-by-N matrix B. On exit, B is overwritten by Q*B or Q**H*B. |
PLASMA_SUCCESS | successful exit |
Definition at line 250 of file zunmqr.c.
References plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_zunmqr_Tile_Async(), and plasma_sequence_t::status.