|
PLASMA
2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
|
This is the group of single complex functions using the advanced synchronous interface.
| int PLASMA_cgebrd_Tile | ( | PLASMA_enum | jobu, |
| PLASMA_enum | jobvt, | ||
| PLASMA_desc * | A, | ||
| float * | D, | ||
| float * | E, | ||
| PLASMA_desc * | U, | ||
| PLASMA_desc * | VT, | ||
| PLASMA_desc * | T | ||
| ) |
PLASMA_cgebrd_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_cgebrd(). 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 cgebrd.c.
References PLASMA_cgebrd_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgelqf_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | T | ||
| ) |
PLASMA_cgelqf_Tile - Computes the tile LQ factorization of a matrix. Tile equivalent of PLASMA_cgelqf(). 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_cgelqs to solve the system of equations. |
| PLASMA_SUCCESS | successful exit |
Definition at line 187 of file cgelqf.c.
References PLASMA_cgelqf_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgelqs_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | T, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgelqs_Tile - Computes a minimum-norm solution using previously computed LQ factorization. Tile equivalent of PLASMA_cgelqs(). 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_cgelqf. |
| [in] | T | Auxiliary factorization data, computed by PLASMA_cgelqf. |
| [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 cgelqs.c.
References PLASMA_cgelqs_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgels_Tile | ( | PLASMA_enum | trans, |
| PLASMA_desc * | A, | ||
| PLASMA_desc * | T, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgels_Tile - Solves overdetermined or underdetermined linear system of equations using the tile QR or the tile LQ factorization. Tile equivalent of PLASMA_cgels(). 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_cgeqrf; if M < N, A is overwritten by details of its LQ factorization as returned by PLASMA_cgelqf. |
| [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 cgels.c.
References PLASMA_cgels_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgemm_Tile | ( | PLASMA_enum | transA, |
| PLASMA_enum | transB, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B, | ||
| PLASMA_Complex32_t | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_cgemm_Tile - Performs matrix multiplication. Tile equivalent of PLASMA_cgemm(). 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 conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is transposed; = PlasmaConjTrans: A is conjfugate transposed. |
| [in] | transB | Specifies whether the matrix B is transposed, not transposed or conjfugate transposed: = PlasmaNoTrans: B is not transposed; = PlasmaTrans: B is transposed; = PlasmaConjTrans: B is conjfugate 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 cgemm.c.
References PLASMA_cgemm_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgeqrf_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | T | ||
| ) |
PLASMA_cgeqrf_Tile - Computes the tile QR factorization of a matrix. Tile equivalent of PLASMA_cgeqrf(). 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_cgeqrs to solve the system of equations. |
| PLASMA_SUCCESS | successful exit |
Definition at line 186 of file cgeqrf.c.
References PLASMA_cgeqrf_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgeqrs_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | T, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgeqrs_Tile - Computes a minimum-norm solution using the tile QR factorization. Tile equivalent of PLASMA_cgeqrf(). 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_cgeqrf. |
| [in] | T | Auxiliary factorization data, computed by PLASMA_cgeqrf. |
| [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 cgeqrs.c.
References PLASMA_cgeqrs_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgesv_incpiv_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | L, | ||
| int * | IPIV, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgesv_incpiv_Tile - Solves a system of linear equations using the tile LU factorization. Tile equivalent of PLASMA_cgetrf_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 cgesv_incpiv.c.
References PLASMA_cgesv_incpiv_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgesv_Tile | ( | PLASMA_desc * | A, |
| int * | IPIV, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgesv_Tile - Solves a system of linear equations using the tile LU factorization. Tile equivalent of PLASMA_cgetrf(). 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 cgesv.c.
References PLASMA_cgesv_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgesvd_Tile | ( | PLASMA_enum | jobu, |
| PLASMA_enum | jobvt, | ||
| PLASMA_desc * | A, | ||
| float * | S, | ||
| PLASMA_desc * | U, | ||
| PLASMA_desc * | VT, | ||
| PLASMA_desc * | T | ||
| ) |
PLASMA_cgesvd_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_cgesvd(). 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 cgesvd.c.
References PLASMA_cgesvd_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgetrf_incpiv_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | L, | ||
| int * | IPIV | ||
| ) |
PLASMA_cgetrf_incpiv_Tile - Computes the tile LU factorization of a matrix. Tile equivalent of PLASMA_cgetrf_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_cgetrs_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 cgetrf_incpiv.c.
References PLASMA_cgetrf_incpiv_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgetrf_Tile | ( | PLASMA_desc * | A, |
| int * | IPIV | ||
| ) |
PLASMA_cgetrf_Tile - Computes the tile LU factorization of a matrix. Tile equivalent of PLASMA_cgetrf(). 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_cgetrs 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 cgetrf.c.
References PLASMA_cgetrf_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgetri_Tile | ( | PLASMA_desc * | A, |
| int * | IPIV | ||
| ) |
PLASMA_cgetri_Tile - Computes the inverse of a matrix using the LU factorization computed by PLASMA_cgetrf. 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_cgetri(). 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_cgetrf. 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_cgetrf. |
| 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 cgetri.c.
References PLASMA_Alloc_Workspace_cgetri_Tile_Async(), PLASMA_cgetri_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(), and plasma_sequence_t::status.

| int PLASMA_cgetrs_incpiv_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | L, | ||
| int * | IPIV, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgetrs_incpiv_Tile - Solves a system of linear equations using previously computed LU factorization. Tile equivalent of PLASMA_cgetrs_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_cgetrf_incpiv. |
| [in] | L | Auxiliary factorization data, related to the tile L factor, computed by PLASMA_cgetrf_incpiv. |
| [in] | IPIV | The pivot indices from PLASMA_cgetrf_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 cgetrs_incpiv.c.
References PLASMA_cgetrs_incpiv_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cgetrs_Tile | ( | PLASMA_enum | trans, |
| PLASMA_desc * | A, | ||
| int * | IPIV, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cgetrs_Tile - Solves a system of linear equations using previously computed LU factorization. Tile equivalent of PLASMA_cgetrs(). 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_cgetrf. |
| [in] | IPIV | The pivot indices from PLASMA_cgetrf. |
| [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 cgetrs.c.
References PLASMA_cgetrs_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cheev_Tile | ( | PLASMA_enum | jobz, |
| PLASMA_enum | uplo, | ||
| PLASMA_desc * | A, | ||
| float * | W, | ||
| PLASMA_desc * | T, | ||
| PLASMA_desc * | Q | ||
| ) |
PLASMA_cheev_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_cheev 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 cheev.c.
References PLASMA_cheev_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_chegst_Tile | ( | PLASMA_enum | itype, |
| PLASMA_enum | uplo, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_chegst_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_CPOTRF. ONLY PlasmaItype == 1 and PlasmaLower supported! Tile equivalent of PLASMA_chegst(). 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_CPOTRF. |
| PLASMA_SUCCESS | successful exit |
| <0 | if -i, the i-th argument had an illegal value |
Definition at line 233 of file chegst.c.
References PLASMA_chegst_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_chegv_Tile | ( | PLASMA_enum | itype, |
| PLASMA_enum | jobz, | ||
| PLASMA_enum | uplo, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B, | ||
| float * | W, | ||
| PLASMA_desc * | T, | ||
| PLASMA_desc * | Q | ||
| ) |
PLASMA_chegv_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_chegv(). 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_chegv 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_chegv 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 chegv.c.
References PLASMA_chegv_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_chemm_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | uplo, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B, | ||
| PLASMA_Complex32_t | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_chemm_Tile - Performs Hermitian matrix multiplication. Tile equivalent of PLASMA_chemm(). 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:
|
| [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 chemm.c.
References PLASMA_chemm_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cher2k_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_enum | trans, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B, | ||
| float | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_cher2k_Tile - Performs hermitian rank k update. Tile equivalent of PLASMA_cher2k(). 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 conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate 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 cher2k.c.
References PLASMA_cher2k_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cherk_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_enum | trans, | ||
| float | alpha, | ||
| PLASMA_desc * | A, | ||
| float | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_cherk_Tile - Performs hermitian rank k update. Tile equivalent of PLASMA_cherk(). 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 conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate 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 cherk.c.
References PLASMA_cherk_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_chetrd_Tile | ( | PLASMA_enum | jobz, |
| PLASMA_enum | uplo, | ||
| PLASMA_desc * | A, | ||
| float * | D, | ||
| float * | E, | ||
| PLASMA_desc * | T, | ||
| PLASMA_desc * | Q | ||
| ) |
PLASMA_chetrd_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_chetrd(). 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 chetrd.c.
References PLASMA_chetrd_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_clacpy_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_clacpy_Tile - Tile equivalent of PLASMA_clacpy(). 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 clacpy.c.
References PLASMA_clacpy_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and PLASMA_SUCCESS.

| float PLASMA_clange_Tile | ( | PLASMA_enum | norm, |
| PLASMA_desc * | A, | ||
| float * | work | ||
| ) |
PLASMA_clange_Tile - Tile equivalent of PLASMA_clange(). 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 | real 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 clange.c.
References PLASMA_clange_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), and plasma_sequence_destroy().

| float PLASMA_clanhe_Tile | ( | PLASMA_enum | norm, |
| PLASMA_enum | uplo, | ||
| PLASMA_desc * | A, | ||
| float * | work | ||
| ) |
PLASMA_clanhe_Tile - Tile equivalent of PLASMA_clanhe(). 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 | real array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise. |
| PLASMA_SUCCESS | successful exit |
Definition at line 195 of file clanhe.c.
References PLASMA_clanhe_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), and plasma_sequence_destroy().

| float PLASMA_clansy_Tile | ( | PLASMA_enum | norm, |
| PLASMA_enum | uplo, | ||
| PLASMA_desc * | A, | ||
| float * | work | ||
| ) |
PLASMA_clansy_Tile - Tile equivalent of PLASMA_clansy(). 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 | real array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise. |
| PLASMA_SUCCESS | successful exit |
Definition at line 195 of file clansy.c.
References PLASMA_clansy_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), and plasma_sequence_destroy().

| int PLASMA_claset_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_Complex32_t | beta, | ||
| PLASMA_desc * | A | ||
| ) |
PLASMA_claset_Tile - Tile equivalent of PLASMA_claset(). 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 claset.c.
References PLASMA_claset_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and PLASMA_SUCCESS.

| int PLASMA_claswp_Tile | ( | PLASMA_desc * | A, |
| int | K1, | ||
| int | K2, | ||
| int * | IPIV, | ||
| int | INCX | ||
| ) |
PLASMA_claswp_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_claswp(). 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_cgetrf. |
| [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_cgetrf. |
| [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 claswp.c.
References PLASMA_claswp_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_claswpc_Tile | ( | PLASMA_desc * | A, |
| int | K1, | ||
| int | K2, | ||
| int * | IPIV, | ||
| int | INCX | ||
| ) |
PLASMA_claswpc_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_claswpc(). 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_cgetrf. |
| [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_cgetrf. |
| [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 claswpc.c.
References PLASMA_claswpc_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_clauum_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_desc * | A | ||
| ) |
PLASMA_clauum_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_clauum(). 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 clauum.c.
References PLASMA_clauum_Tile_Async(), plasma_context_self(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cplghe_Tile | ( | float | bump, |
| PLASMA_desc * | A, | ||
| unsigned long long int | seed | ||
| ) |
PLASMA_cplghe_Tile - Generate a random hermitian matrix by tiles. Tile equivalent of PLASMA_cplghe(). 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 cplghe.c.
References plasma_context_self(), PLASMA_cplghe_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cplgsy_Tile | ( | PLASMA_Complex32_t | bump, |
| PLASMA_desc * | A, | ||
| unsigned long long int | seed | ||
| ) |
PLASMA_cplgsy_Tile - Generate a random hermitian matrix by tiles. Tile equivalent of PLASMA_cplgsy(). 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 cplgsy.c.
References plasma_context_self(), PLASMA_cplgsy_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cplrnt_Tile | ( | PLASMA_desc * | A, |
| unsigned long long int | seed | ||
| ) |
PLASMA_cplrnt_Tile - Generate a random matrix by tiles. Tile equivalent of PLASMA_cplrnt(). 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 cplrnt.c.
References plasma_context_self(), PLASMA_cplrnt_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cposv_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cposv_Tile - Solves a symmetric positive definite or Hermitian positive definite system of linear equations using the Cholesky factorization. Tile equivalent of PLASMA_cposv(). 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 cposv.c.
References plasma_context_self(), PLASMA_cposv_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cpotrf_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_desc * | A | ||
| ) |
PLASMA_cpotrf_Tile - Computes the Cholesky factorization of a symmetric positive definite or Hermitian positive definite matrix. Tile equivalent of PLASMA_cpotrf(). 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 cpotrf.c.
References plasma_context_self(), PLASMA_cpotrf_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cpotri_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_desc * | A | ||
| ) |
PLASMA_cpotri_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_cpotrf. Tile equivalent of PLASMA_cpotri(). 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_cpotrf. 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 cpotri.c.
References plasma_context_self(), PLASMA_cpotri_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cpotrs_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cpotrs_Tile - Solves a system of linear equations using previously computed Cholesky factorization. Tile equivalent of PLASMA_cpotrs(). 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_cpotrf. |
| [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 cpotrs.c.
References plasma_context_self(), PLASMA_cpotrs_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_csymm_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | uplo, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B, | ||
| PLASMA_Complex32_t | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_csymm_Tile - Performs symmetric matrix multiplication. Tile equivalent of PLASMA_csymm(). 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:
|
| [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 csymm.c.
References plasma_context_self(), PLASMA_csymm_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_csyr2k_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_enum | trans, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B, | ||
| PLASMA_Complex32_t | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_csyr2k_Tile - Performs symmetric rank k update. Tile equivalent of PLASMA_csyr2k(). 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 conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is conjfugate 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 csyr2k.c.
References plasma_context_self(), PLASMA_csyr2k_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_csyrk_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_enum | trans, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_Complex32_t | beta, | ||
| PLASMA_desc * | C | ||
| ) |
PLASMA_csyrk_Tile - Performs rank k update. Tile equivalent of PLASMA_csyrk(). 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 conjfugate 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 csyrk.c.
References plasma_context_self(), PLASMA_csyrk_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_ctrmm_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | uplo, | ||
| PLASMA_enum | transA, | ||
| PLASMA_enum | diag, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_ctrmm_Tile - Computes triangular solve. Tile equivalent of PLASMA_ctrmm(). 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 conjfugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate 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 ctrmm.c.
References plasma_context_self(), PLASMA_ctrmm_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_ctrsm_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | uplo, | ||
| PLASMA_enum | transA, | ||
| PLASMA_enum | diag, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_ctrsm_Tile - Computes triangular solve. Tile equivalent of PLASMA_ctrsm(). 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 conjfugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate 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 ctrsm.c.
References plasma_context_self(), PLASMA_ctrsm_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_ctrsmpl_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | L, | ||
| int * | IPIV, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_ctrsmpl_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_cgetrf_incpiv. |
| [in] | L | Auxiliary factorization data, related to the tile L factor, computed by PLASMA_cgetrf_incpiv. |
| [in] | IPIV | The pivot indices from PLASMA_cgetrf_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 ctrsmpl.c.
References plasma_context_self(), PLASMA_ctrsmpl_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_ctrsmrv_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | uplo, | ||
| PLASMA_enum | transA, | ||
| PLASMA_enum | diag, | ||
| PLASMA_Complex32_t | alpha, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_ctrsmrv_Tile - Computes triangular solve. Tile equivalent of PLASMA_ctrsmrv(). 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 conjfugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate 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 ctrsmrv.c.
References plasma_context_self(), PLASMA_ctrsmrv_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_ctrtri_Tile | ( | PLASMA_enum | uplo, |
| PLASMA_enum | diag, | ||
| PLASMA_desc * | A | ||
| ) |
PLASMA_ctrtri_Tile - Computes the inverse of a complex upper or lower triangular matrix A. Tile equivalent of PLASMA_ctrtri(). 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 ctrtri.c.
References plasma_context_self(), PLASMA_ctrtri_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

| int PLASMA_cunglq_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | T, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cunglq_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_cgelqf. 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_cgelqf. |
| [in] | T | Auxiliary factorization data, computed by PLASMA_cgelqf. |
| [out] | B | On exit, the M-by-N matrix Q. |
| PLASMA_SUCCESS | successful exit |
Definition at line 202 of file cunglq.c.
References plasma_context_self(), PLASMA_cunglq_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cungqr_Tile | ( | PLASMA_desc * | A, |
| PLASMA_desc * | T, | ||
| PLASMA_desc * | Q | ||
| ) |
PLASMA_cungqr_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_cgeqrf. 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_cgeqrf. |
| [in] | T | Auxiliary factorization data, computed by PLASMA_cgeqrf. |
| [out] | Q | On exit, the M-by-N matrix Q. |
| PLASMA_SUCCESS | successful exit |
Definition at line 200 of file cungqr.c.
References plasma_context_self(), PLASMA_cungqr_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cunmlq_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | trans, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | T, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cunmlq_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_cgelqf_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: conjfugate 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_cgelqf. |
| [in] | T | Auxiliary factorization data, computed by PLASMA_cgelqf. |
| [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 cunmlq.c.
References plasma_context_self(), PLASMA_cunmlq_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.


| int PLASMA_cunmqr_Tile | ( | PLASMA_enum | side, |
| PLASMA_enum | trans, | ||
| PLASMA_desc * | A, | ||
| PLASMA_desc * | T, | ||
| PLASMA_desc * | B | ||
| ) |
PLASMA_cunmqr_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_cgeqrf_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: conjfugate 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_cgeqrf. |
| [in] | T | Auxiliary factorization data, computed by PLASMA_cgeqrf. |
| [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 cunmqr.c.
References plasma_context_self(), PLASMA_cunmqr_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), and plasma_sequence_t::status.

