PLASMA  2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
Advanced Interface: Synchronous - Single Complex

Functions

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)
int PLASMA_cgelqf_Tile (PLASMA_desc *A, PLASMA_desc *T)
int PLASMA_cgelqs_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_cgels_Tile (PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
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)
int PLASMA_cgeqrf_Tile (PLASMA_desc *A, PLASMA_desc *T)
int PLASMA_cgeqrs_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_cgesv_Tile (PLASMA_desc *A, int *IPIV, PLASMA_desc *B)
int PLASMA_cgesv_incpiv_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV, PLASMA_desc *B)
int PLASMA_cgesvd_Tile (PLASMA_enum jobu, PLASMA_enum jobvt, PLASMA_desc *A, float *S, PLASMA_desc *U, PLASMA_desc *VT, PLASMA_desc *T)
int PLASMA_cgetrf_Tile (PLASMA_desc *A, int *IPIV)
int PLASMA_cgetrf_incpiv_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV)
int PLASMA_cgetri_Tile (PLASMA_desc *A, int *IPIV)
int PLASMA_cgetrs_Tile (PLASMA_enum trans, PLASMA_desc *A, int *IPIV, PLASMA_desc *B)
int PLASMA_cgetrs_incpiv_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV, PLASMA_desc *B)
int PLASMA_cheev_Tile (PLASMA_enum jobz, PLASMA_enum uplo, PLASMA_desc *A, float *W, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_chegst_Tile (PLASMA_enum itype, PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
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)
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)
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)
int PLASMA_cherk_Tile (PLASMA_enum uplo, PLASMA_enum trans, float alpha, PLASMA_desc *A, float beta, PLASMA_desc *C)
int PLASMA_chetrd_Tile (PLASMA_enum jobz, PLASMA_enum uplo, PLASMA_desc *A, float *D, float *E, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_clacpy_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
float PLASMA_clange_Tile (PLASMA_enum norm, PLASMA_desc *A, float *work)
float PLASMA_clanhe_Tile (PLASMA_enum norm, PLASMA_enum uplo, PLASMA_desc *A, float *work)
float PLASMA_clansy_Tile (PLASMA_enum norm, PLASMA_enum uplo, PLASMA_desc *A, float *work)
int PLASMA_claset_Tile (PLASMA_enum uplo, PLASMA_Complex32_t alpha, PLASMA_Complex32_t beta, PLASMA_desc *A)
int PLASMA_claswp_Tile (PLASMA_desc *A, int K1, int K2, int *IPIV, int INCX)
int PLASMA_claswpc_Tile (PLASMA_desc *A, int K1, int K2, int *IPIV, int INCX)
int PLASMA_clauum_Tile (PLASMA_enum uplo, PLASMA_desc *A)
int PLASMA_cplghe_Tile (float bump, PLASMA_desc *A, unsigned long long int seed)
int PLASMA_cplgsy_Tile (PLASMA_Complex32_t bump, PLASMA_desc *A, unsigned long long int seed)
int PLASMA_cplrnt_Tile (PLASMA_desc *A, unsigned long long int seed)
int PLASMA_cposv_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_cpotrf_Tile (PLASMA_enum uplo, PLASMA_desc *A)
int PLASMA_cpotri_Tile (PLASMA_enum uplo, PLASMA_desc *A)
int PLASMA_cpotrs_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
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)
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)
int PLASMA_csyrk_Tile (PLASMA_enum uplo, PLASMA_enum trans, PLASMA_Complex32_t alpha, PLASMA_desc *A, PLASMA_Complex32_t beta, PLASMA_desc *C)
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)
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)
int PLASMA_ctrsmpl_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV, PLASMA_desc *B)
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)
int PLASMA_ctrtri_Tile (PLASMA_enum uplo, PLASMA_enum diag, PLASMA_desc *A)
int PLASMA_cunglq_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_cungqr_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_cunmlq_Tile (PLASMA_enum side, PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_cunmqr_Tile (PLASMA_enum side, PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)

Detailed Description

This is the group of single complex functions using the advanced synchronous interface.


Function Documentation

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.

Parameters:
[in]jobuSpecifies 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]jobvtSpecifies 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]AOn 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]SThe 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]VTIf 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]TOn exit, auxiliary factorization data.
Returns:
PLASMA_SUCCESS successful exit
See also:
PLASMA_cgebrd
PLASMA_cgebrd_Tile_Async
PLASMA_cgebrd_Tile
PLASMA_dgebrd_Tile
PLASMA_sgebrd_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgebrd_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgebrd_Tile_Async(jobu, jobvt, A, D, E, U, VT, T, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn 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]TOn exit, auxiliary factorization data, required by PLASMA_cgelqs to solve the system of equations.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgelqf
PLASMA_cgelqf_Tile_Async
PLASMA_cgelqf_Tile
PLASMA_dgelqf_Tile
PLASMA_sgelqf_Tile
PLASMA_cgelqs_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgelqf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgelqf_Tile_Async(A, T, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]ADetails of the LQ factorization of the original matrix A as returned by PLASMA_cgelqf.
[in]TAuxiliary factorization data, computed by PLASMA_cgelqf.
[in,out]BOn entry, the M-by-NRHS right hand side matrix B. On exit, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgelqs
PLASMA_cgelqs_Tile_Async
PLASMA_cgelqs_Tile
PLASMA_dgelqs_Tile
PLASMA_sgelqs_Tile
PLASMA_cgelqf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgelqs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgelqs_Tile_Async(A, T, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]transIntended usage: = PlasmaNoTrans: the linear system involves A; = PlasmaConjTrans: the linear system involves A**H. Currently only PlasmaNoTrans is supported.
[in,out]AOn 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]TOn exit, auxiliary factorization data.
[in,out]BOn 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;
Returns:
PLASMA_SUCCESS successful exit
See also:
PLASMA_cgels
PLASMA_cgels_Tile_Async
PLASMA_cgels_Tile
PLASMA_dgels_Tile
PLASMA_sgels_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgels_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgels_Tile_Async(trans, A, T, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]transASpecifies 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]transBSpecifies 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]alphaalpha specifies the scalar alpha
[in]AA is a LDA-by-ka matrix, where ka is K when transA = PlasmaNoTrans, and is M otherwise.
[in]BB is a LDB-by-kb matrix, where kb is N when transB = PlasmaNoTrans, and is K otherwise.
[in]betabeta specifies the scalar beta
[in,out]CC 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 )
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgemm
PLASMA_cgemm_Tile_Async
PLASMA_cgemm_Tile
PLASMA_dgemm_Tile
PLASMA_sgemm_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgemm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgemm_Tile_Async(transA, transB, alpha, A, B, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn 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]TOn exit, auxiliary factorization data, required by PLASMA_cgeqrs to solve the system of equations.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgeqrf
PLASMA_cgeqrf_Tile_Async
PLASMA_cgeqrf_Tile
PLASMA_dgeqrf_Tile
PLASMA_sgeqrf_Tile
PLASMA_cgeqrs_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgeqrf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgeqrf_Tile_Async(A, T, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]ADetails of the QR factorization of the original matrix A as returned by PLASMA_cgeqrf.
[in]TAuxiliary factorization data, computed by PLASMA_cgeqrf.
[in,out]BOn entry, the m-by-nrhs right hand side matrix B. On exit, the n-by-nrhs solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgeqrs
PLASMA_cgeqrs_Tile_Async
PLASMA_cgeqrs_Tile
PLASMA_dgeqrs_Tile
PLASMA_sgeqrs_Tile
PLASMA_cgeqrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgeqrs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgeqrs_Tile_Async(A, T, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn 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]LOn exit, auxiliary factorization data, related to the tile L factor, necessary to solve the system of equations.
[out]IPIVOn exit, the pivot indices that define the permutations (not equivalent to LAPACK).
[in,out]BOn 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if 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.
See also:
PLASMA_cgesv_incpiv
PLASMA_cgesv_incpiv_Tile_Async
PLASMA_cgesv_incpiv_Tile
PLASMA_dgesv_incpiv_Tile
PLASMA_sgesv_incpiv_Tile
PLASMA_ccgesv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgesv_incpiv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgesv_incpiv_Tile_Async(A, L, IPIV, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn entry, the N-by-N coefficient matrix A. On exit, the tile L and U factors from the factorization.
[out]IPIVOn exit, the pivot indices that define the permutations.
[in,out]BOn 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if 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.
See also:
PLASMA_cgesv
PLASMA_cgesv_Tile_Async
PLASMA_cgesv_Tile
PLASMA_dgesv_Tile
PLASMA_sgesv_Tile
PLASMA_ccgesv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgesv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgesv_Tile_Async(A, IPIV, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]jobuSpecifies 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]jobvtSpecifies 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]AOn 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]SThe 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]VTIf 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]TOn exit, auxiliary factorization data.
Returns:
PLASMA_SUCCESS successful exit
See also:
PLASMA_cgesvd
PLASMA_cgesvd_Tile_Async
PLASMA_cgesvd_Tile
PLASMA_dgesvd_Tile
PLASMA_sgesvd_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgesvd_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgesvd_Tile_Async(jobu, jobvt, A, S, U, VT, T, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn entry, the M-by-N matrix to be factored. On exit, the tile factors L and U from the factorization.
[out]LOn exit, auxiliary factorization data, related to the tile L factor, required by PLASMA_cgetrs_incpiv to solve the system of equations.
[out]IPIVThe pivot indices that define the permutations (not equivalent to LAPACK).
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if 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.
See also:
PLASMA_cgetrf_incpiv
PLASMA_cgetrf_incpiv_Tile_Async
PLASMA_cgetrf_incpiv_Tile
PLASMA_dgetrf_incpiv_Tile
PLASMA_sgetrf_incpiv_Tile
PLASMA_cgetrs_incpiv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgetrf_incpiv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgetrf_incpiv_Tile_Async(A, L, IPIV, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn entry, the M-by-N matrix to be factored. On exit, the tile factors L and U from the factorization.
[out]LOn exit, auxiliary factorization data, related to the tile L factor, required by PLASMA_cgetrs to solve the system of equations.
[out]IPIVThe pivot indices that define the permutations (not equivalent to LAPACK).
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if 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.
See also:
PLASMA_cgetrf
PLASMA_cgetrf_Tile_Async
PLASMA_cgetrf_Tile
PLASMA_dgetrf_Tile
PLASMA_sgetrf_Tile
PLASMA_cgetrs_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgetrf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgetrf_Tile_Async(A, IPIV, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in,out]AOn 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]IPIVThe pivot indices that define the permutations as returned by PLASMA_cgetrf.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if i, the (i,i) element of the factor U is exactly zero; The matrix is singular and its inverse could not be computed.
See also:
PLASMA_cgetri
PLASMA_cgetri_Tile_Async
PLASMA_cgetri_Tile
PLASMA_dgetri_Tile
PLASMA_sgetri_Tile
PLASMA_cgetrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
PLASMA_desc descW;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgetri_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
/* Allocate workspace */
PLASMA_cgetri_Tile_Async(A, IPIV, &descW, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]AThe tile factors L and U from the factorization, computed by PLASMA_cgetrf_incpiv.
[in]LAuxiliary factorization data, related to the tile L factor, computed by PLASMA_cgetrf_incpiv.
[in]IPIVThe pivot indices from PLASMA_cgetrf_incpiv (not equivalent to LAPACK).
[in,out]BOn entry, the N-by-NRHS matrix of right hand side matrix B. On exit, the solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgetrs_incpiv
PLASMA_cgetrs_incpiv_Tile_Async
PLASMA_cgetrs_incpiv_Tile
PLASMA_dgetrs_incpiv_Tile
PLASMA_sgetrs_incpiv_Tile
PLASMA_cgetrf_incpiv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgetrs_incpiv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgetrs_incpiv_Tile_Async(A, L, IPIV, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]transIntended 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]AThe tile factors L and U from the factorization, computed by PLASMA_cgetrf.
[in]IPIVThe pivot indices from PLASMA_cgetrf.
[in,out]BOn entry, the N-by-NRHS matrix of right hand side matrix B. On exit, the solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cgetrs
PLASMA_cgetrs_Tile_Async
PLASMA_cgetrs_Tile
PLASMA_dgetrs_Tile
PLASMA_sgetrs_Tile
PLASMA_cgetrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cgetrs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cgetrs_Tile_Async(trans, A, IPIV, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]jobzIntended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. PlasmaVec is NOT supported.
[in]uploSpecifies 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]AOn 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]WOn exit, if info = 0, the eigenvalues.
[in,out]TOn entry, descriptor as return by PLASMA_Alloc_Workspace_cheev On exit, contains auxiliary factorization data.
[out]QOn exit, if jobz = PlasmaVec and info = 0, the eigenvectors.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value
>0if INFO = i, the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
See also:
PLASMA_cheev_Tile
PLASMA_cheev_Tile_Async
PLASMA_cheev
PLASMA_dsyev
PLASMA_ssyev

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cheev_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cheev_Tile_Async(jobz, uplo, A, W, T, Q, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]PlasmaItypeIntended 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]uploSpecifies 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]AOn 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]BOn entry, the triangular factor from the Cholesky factorization of B, as returned by PLASMA_CPOTRF.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value
See also:
PLASMA_chegst
PLASMA_chegst_Tile_Async
PLASMA_chegst_Tile
PLASMA_dsygst_Tile
PLASMA_ssygst_Tile
PLASMA_chegst_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_chegst_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_chegst_Tile_Async(itype, uplo, A, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]PlasmaItypeIntended usage: = 1: A*x=(lambda)*B*x = 2: A*Bx=(lambda)*x = 3: B*A*x=(lambda)*x
[in]jobzIntended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. Currently only PlasmaVec is NOT supported.
[in]uploSpecifies 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]NThe order of the matrix A. N >= 0.
[in,out]AOn 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]LDAThe leading dimension of the array A. LDA >= max(1,N).
[in,out]BOn 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]LDBThe leading dimension of the array B. LDA >= max(1,N).
[in,out]TOn entry, descriptor as return by PLASMA_Alloc_Workspace_chegv On exit, contains auxiliary factorization data.
[out]WOn exit, if info = 0, the eigenvalues.
[out]QOn exit, if jobz = PlasmaVec and info = 0, the eigenvectors.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value
<=Nif INFO = i, plasma_chegv failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
>Nif 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.
See also:
PLASMA_chegv
PLASMA_chegv_Tile_Async
PLASMA_chegv_Tile
PLASMA_dsygv_Tile
PLASMA_ssygv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_chegv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_chegv_Tile_Async(itype, jobz, uplo, A, B, W, T, Q, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]sideSpecifies whether the hermitian matrix A appears on the left or right in the operation as follows: = PlasmaLeft:

\[ C = \alpha \times A \times B + \beta \times C \]

= PlasmaRight:

\[ C = \alpha \times B \times A + \beta \times C \]

[in]uploSpecifies 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]alphaSpecifies the scalar alpha.
[in]AA 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]BB is a LDB-by-N matrix, where the leading M-by-N part of the array B must contain the matrix B.
[in]betaSpecifies the scalar beta.
[in,out]CC is a LDC-by-N matrix. On exit, the array is overwritten by the M by N updated matrix.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_chemm
PLASMA_chemm_Tile_Async
PLASMA_chemm_Tile
PLASMA_dhemm_Tile
PLASMA_shemm_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_chemm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_chemm_Tile_Async(side, uplo, alpha, A, B, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored.
[in]transSpecifies whether the matrix A is transposed or conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate transposed.
[in]alphaalpha specifies the scalar alpha.
[in]AA is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise.
[in]BB is a LDB-by-kb matrix, where kb is K when trans = PlasmaNoTrans, and is N otherwise.
[in]betabeta specifies the scalar beta
[in,out]CC 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cher2k_Tile
PLASMA_cher2k
PLASMA_dher2k
PLASMA_sher2k

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cher2k_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cher2k_Tile_Async(uplo, trans, alpha, A, B, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored.
[in]transSpecifies whether the matrix A is transposed or conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaConjTrans: A is conjfugate transposed.
[in]alphaalpha specifies the scalar alpha.
[in]AA is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise.
[in]betabeta specifies the scalar beta
[in,out]CC 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cherk_Tile
PLASMA_cherk
PLASMA_dherk
PLASMA_sherk

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cherk_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cherk_Tile_Async(uplo, trans, alpha, A, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]jobzIntended usage: = PlasmaNoVec: computes eigenvalues only; = PlasmaVec: computes eigenvalues and eigenvectors. PlasmaVec is NOT supported.
[in]uploSpecifies 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]AOn 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]DOn exit, the diagonal elements of the tridiagonal matrix: D(i) = A(i,i).
[out]EOn 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]TOn exit, auxiliary factorization data.
[out]QOn exit, if jobz = PlasmaVec and info = 0, the eigenvectors.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value
>0if INFO = i, the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
See also:
PLASMA_chetrd
PLASMA_chetrd_Tile_Async
PLASMA_chetrd_Tile
PLASMA_dsytrd_Tile
PLASMA_ssytrd_Tile
PLASMA_cheev_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_chetrd_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_chetrd_Tile_Async(jobz, uplo, A, D, E, T, Q, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]uploSpecifies 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]AThe M-by-N matrix A. If uplo = PlasmaUpper, only the upper trapezium is accessed; if UPLO = PlasmaLower, only the lower trapezium is accessed.
[out]BThe M-by-N matrix B. On exit, B = A in the locations specified by UPLO.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_clacpy
PLASMA_clacpy_Tile_Async
PLASMA_clacpy_Tile
PLASMA_dlacpy_Tile
PLASMA_slacpy_Tile

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.

{
PLASMA_sequence *sequence = NULL;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_clacpy_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_clacpy_Tile_Async(uplo, A, B, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
}

Here is the call graph for this function:

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.

Parameters:
[in]norm= PlasmaMaxNorm: Max norm = PlasmaOneNorm: One norm = PlasmaInfNorm: Infinity norm = PlasmaFrobeniusNorm: Frobenius norm
[in]AOn 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]workreal array of dimension (MAX(1,LWORK)), where LWORK >= M when NORM = PlasmaInfNorm; otherwise, WORK is not referenced.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_clange
PLASMA_clange_Tile_Async
PLASMA_clange_Tile
PLASMA_dlange_Tile
PLASMA_slange_Tile

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().

{
PLASMA_sequence *sequence = NULL;
float value;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_clange_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_clange_Tile_Async(norm, A, work, &value, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
return value;
}

Here is the call graph for this function:

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.

Parameters:
[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]AOn 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]workreal array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_clanhe
PLASMA_clanhe_Tile_Async
PLASMA_clanhe_Tile
PLASMA_dlanhe_Tile
PLASMA_slanhe_Tile

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().

{
PLASMA_sequence *sequence = NULL;
float value;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_clanhe_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_clanhe_Tile_Async(norm, uplo, A, work, &value, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
return value;
}

Here is the call graph for this function:

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.

Parameters:
[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]AOn 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]workreal array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_clansy
PLASMA_clansy_Tile_Async
PLASMA_clansy_Tile
PLASMA_dlansy_Tile
PLASMA_slansy_Tile

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().

{
PLASMA_sequence *sequence = NULL;
float value;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_clansy_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_clansy_Tile_Async(norm, uplo, A, work, &value, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
return value;
}

Here is the call graph for this function:

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.

Parameters:
[in]uploSpecifies 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]AOn 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)
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_claset
PLASMA_claset_Tile_Async
PLASMA_claset_Tile
PLASMA_dlaset_Tile
PLASMA_slaset_Tile

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.

{
PLASMA_sequence *sequence = NULL;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_claset_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_claset_Tile_Async(uplo, alpha, beta, A, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
}

Here is the call graph for this function:

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.

Parameters:
[in]AThe tile factors L and U from the factorization, computed by PLASMA_cgetrf.
[in]K1The first element of IPIV for which a row interchange will be done.
[in]K2The last element of IPIV for which a row interchange will be done.
[in]IPIVThe pivot indices from PLASMA_cgetrf.
[in]INCXThe increment between successive values of IPIV. If IPIV is negative, the pivots are applied in reverse order.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_claswp
PLASMA_claswp_Tile_Async
PLASMA_claswp_Tile
PLASMA_dlaswp_Tile
PLASMA_slaswp_Tile
PLASMA_cgetrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_claswp_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_claswp_Tile_Async(A, K1, K2, IPIV, INCX, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]AThe tile factors L and U from the factorization, computed by PLASMA_cgetrf.
[in]K1The first element of IPIV for which a row interchange will be done.
[in]K2The last element of IPIV for which a row interchange will be done.
[in]IPIVThe pivot indices from PLASMA_cgetrf.
[in]INCXThe increment between successive values of IPIV. If IPIV is negative, the pivots are applied in reverse order.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_claswpc
PLASMA_claswpc_Tile_Async
PLASMA_claswpc_Tile
PLASMA_dlaswpc_Tile
PLASMA_slaswpc_Tile
PLASMA_cgetrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_claswpc_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_claswpc_Tile_Async(A, K1, K2, IPIV, INCX, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored.
[in]AOn 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_clauum
PLASMA_clauum_Tile_Async
PLASMA_clauum_Tile
PLASMA_dlauum_Tile
PLASMA_slauum_Tile
PLASMA_cpotri_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_clauum_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_clauum_Tile_Async(uplo, A, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]bumpThe value to add to the diagonal to be sure to have a positive definite matrix.
[in]AOn exit, The random hermitian matrix A generated.
[in]seedThe seed used in the random generation.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cplghe
PLASMA_cplghe_Tile_Async
PLASMA_cplghe_Tile
PLASMA_dplghe_Tile
PLASMA_splghe_Tile
PLASMA_cplrnt_Tile
PLASMA_cplgsy_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cplghe_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cplghe_Tile_Async( bump, A, seed, sequence, &request );
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]bumpThe value to add to the diagonal to be sure to have a positive definite matrix.
[in]AOn exit, The random hermitian matrix A generated.
[in]seedThe seed used in the random generation.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cplgsy
PLASMA_cplgsy_Tile_Async
PLASMA_cplgsy_Tile
PLASMA_dplgsy_Tile
PLASMA_splgsy_Tile
PLASMA_cplrnt_Tile
PLASMA_cplgsy_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cplgsy_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cplgsy_Tile_Async( bump, A, seed, sequence, &request );
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]AOn exit, The random matrix A generated.
[in]seedThe seed used in the random generation.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cplrnt
PLASMA_cplrnt_Tile_Async
PLASMA_cplrnt_Tile
PLASMA_dplrnt_Tile
PLASMA_splrnt_Tile
PLASMA_cplghe_Tile
PLASMA_cplgsy_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cplrnt_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cplrnt_Tile_Async( A, seed, sequence, &request );
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]uploSpecifies 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]AOn 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]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if i, the leading minor of order i of A is not positive definite, so the factorization could not be completed, and the solution has not been computed.
See also:
PLASMA_cposv
PLASMA_cposv_Tile_Async
PLASMA_cposv_Tile
PLASMA_dposv_Tile
PLASMA_sposv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cposv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cposv_Tile_Async(uplo, A, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored.
[in]AOn 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if i, the leading minor of order i of A is not positive definite, so the factorization could not be completed, and the solution has not been computed.
See also:
PLASMA_cpotrf
PLASMA_cpotrf_Tile_Async
PLASMA_cpotrf_Tile
PLASMA_dpotrf_Tile
PLASMA_spotrf_Tile
PLASMA_cpotrs_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cpotrf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cpotrf_Tile_Async(uplo, A, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored.
[in]AOn 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if i, the leading minor of order i of A is not positive definite, so the factorization could not be completed, and the solution has not been computed.
See also:
PLASMA_cpotri
PLASMA_cpotri_Tile_Async
PLASMA_cpotri_Tile
PLASMA_dpotri_Tile
PLASMA_spotri_Tile
PLASMA_cpotrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cpotri_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cpotri_Tile_Async(uplo, A, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of A is stored; = PlasmaLower: Lower triangle of A is stored.
[in]AThe triangular factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H, computed by PLASMA_cpotrf.
[in,out]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cpotrs
PLASMA_cpotrs_Tile_Async
PLASMA_cpotrs_Tile
PLASMA_dpotrs_Tile
PLASMA_spotrs_Tile
PLASMA_cpotrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cpotrs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cpotrs_Tile_Async(uplo, A, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]sideSpecifies whether the symmetric matrix A appears on the left or right in the operation as follows: = PlasmaLeft:

\[ C = \alpha \times A \times B + \beta \times C \]

= PlasmaRight:

\[ C = \alpha \times B \times A + \beta \times C \]

[in]uploSpecifies 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]alphaSpecifies the scalar alpha.
[in]AA 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]BB is a LDB-by-N matrix, where the leading M-by-N part of the array B must contain the matrix B.
[in]betaSpecifies the scalar beta.
[in,out]CC is a LDC-by-N matrix. On exit, the array is overwritten by the M by N updated matrix.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_csymm
PLASMA_csymm_Tile_Async
PLASMA_csymm_Tile
PLASMA_dsymm_Tile
PLASMA_ssymm_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_csymm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_csymm_Tile_Async(side, uplo, alpha, A, B, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored.
[in]transSpecifies whether the matrix A is transposed or conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is conjfugate transposed.
[in]alphaalpha specifies the scalar alpha.
[in]AA is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise.
[in]BB is a LDB-by-kb matrix, where kb is K when trans = PlasmaNoTrans, and is N otherwise.
[in]betabeta specifies the scalar beta
[in,out]CC 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_csyr2k_Tile
PLASMA_csyr2k
PLASMA_dsyr2k
PLASMA_ssyr2k

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_csyr2k_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_csyr2k_Tile_Async(uplo, trans, alpha, A, B, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]uplo= PlasmaUpper: Upper triangle of C is stored; = PlasmaLower: Lower triangle of C is stored.
[in]transSpecifies whether the matrix A is transposed or conjfugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is transposed.
[in]alphaalpha specifies the scalar alpha.
[in]AA is a LDA-by-ka matrix, where ka is K when trans = PlasmaNoTrans, and is N otherwise.
[in]betabeta specifies the scalar beta
[in,out]CC 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_csyrk_Tile
PLASMA_csyrk
PLASMA_dsyrk
PLASMA_ssyrk

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_csyrk_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_csyrk_Tile_Async(uplo, trans, alpha, A, beta, C, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]sideSpecifies whether A appears on the left or on the right of X: = PlasmaLeft: A*X = B = PlasmaRight: X*A = B
[in]uploSpecifies 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]transASpecifies 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]diagSpecifies whether or not A is unit triangular: = PlasmaNonUnit: A is non unit; = PlasmaUnit: A us unit.
[in]alphaalpha specifies the scalar alpha.
[in]AThe 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]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_ctrmm
PLASMA_ctrmm_Tile_Async
PLASMA_ctrmm_Tile
PLASMA_dtrmm_Tile
PLASMA_strmm_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_ctrmm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_ctrmm_Tile_Async(side, uplo, transA, diag, alpha, A, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]sideSpecifies whether A appears on the left or on the right of X: = PlasmaLeft: A*X = B = PlasmaRight: X*A = B
[in]uploSpecifies 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]transASpecifies 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]diagSpecifies whether or not A is unit triangular: = PlasmaNonUnit: A is non unit; = PlasmaUnit: A us unit.
[in]alphaalpha specifies the scalar alpha.
[in]AThe 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]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_ctrsm
PLASMA_ctrsm_Tile_Async
PLASMA_ctrsm_Tile
PLASMA_dtrsm_Tile
PLASMA_strsm_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_ctrsm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_ctrsm_Tile_Async(side, uplo, transA, diag, alpha, A, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]AThe tile factor L from the factorization, computed by PLASMA_cgetrf_incpiv.
[in]LAuxiliary factorization data, related to the tile L factor, computed by PLASMA_cgetrf_incpiv.
[in]IPIVThe pivot indices from PLASMA_cgetrf_incpiv (not equivalent to LAPACK).
[in,out]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_ctrsmpl
PLASMA_ctrsmpl_Tile_Async
PLASMA_ctrsmpl_Tile
PLASMA_dtrsmpl_Tile
PLASMA_strsmpl_Tile
PLASMA_cgetrf_incpiv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_ctrsmpl_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_ctrsmpl_Tile_Async(A, L, IPIV, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]sideSpecifies whether A appears on the left or on the right of X: = PlasmaLeft: A*X = B = PlasmaRight: X*A = B
[in]uploSpecifies 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]transASpecifies 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]diagSpecifies whether or not A is unit triangular: = PlasmaNonUnit: A is non unit; = PlasmaUnit: A us unit.
[in]alphaalpha specifies the scalar alpha.
[in]AThe 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]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_ctrsmrv
PLASMA_ctrsmrv_Tile_Async
PLASMA_ctrsmrv_Tile
PLASMA_dtrsmrv_Tile
PLASMA_strsmrv_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_ctrsmrv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_ctrsmrv_Tile_Async(side, uplo, transA, diag, alpha, A, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[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]AOn 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.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
>0if i, A(i,i) is exactly zero. The triangular matrix is singular and its inverse can not be computed.
See also:
PLASMA_ctrtri
PLASMA_ctrtri_Tile_Async
PLASMA_ctrtri_Tile
PLASMA_dtrtri_Tile
PLASMA_strtri_Tile
PLASMA_cpotri_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_ctrtri_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_ctrtri_Tile_Async(uplo, diag, A, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

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.

Parameters:
[in]ADetails of the LQ factorization of the original matrix A as returned by PLASMA_cgelqf.
[in]TAuxiliary factorization data, computed by PLASMA_cgelqf.
[out]BOn exit, the M-by-N matrix Q.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cunglq
PLASMA_cunglq_Tile_Async
PLASMA_cunglq_Tile
PLASMA_dorglq_Tile
PLASMA_sorglq_Tile
PLASMA_cgelqf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cunglq_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cunglq_Tile_Async(A, T, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]ADetails of the QR factorization of the original matrix A as returned by PLASMA_cgeqrf.
[in]TAuxiliary factorization data, computed by PLASMA_cgeqrf.
[out]QOn exit, the M-by-N matrix Q.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cungqr
PLASMA_cungqr_Tile_Async
PLASMA_cungqr_Tile
PLASMA_dorgqr_Tile
PLASMA_sorgqr_Tile
PLASMA_cgeqrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cungqr_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cungqr_Tile_Async(A, T, Q, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]sideIntended 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]transIntended usage: = PlasmaNoTrans: no transpose, apply Q; = PlasmaConjTrans: conjfugate transpose, apply Q**H. Currently only PlasmaConjTrans is supported.
[in]ADetails of the LQ factorization of the original matrix A as returned by PLASMA_cgelqf.
[in]TAuxiliary factorization data, computed by PLASMA_cgelqf.
[in,out]BOn entry, the M-by-N matrix B. On exit, B is overwritten by Q*B or Q**H*B.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cunmlq
PLASMA_cunmlq_Tile_Async
PLASMA_cunmlq_Tile
PLASMA_dormlq_Tile
PLASMA_sormlq_Tile
PLASMA_cgelqf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cunmlq_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cunmlq_Tile_Async(side, trans, A, T, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Parameters:
[in]sideIntended 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]transIntended usage: = PlasmaNoTrans: no transpose, apply Q; = PlasmaConjTrans: conjfugate transpose, apply Q**H. Currently only PlasmaConjTrans is supported.
[in]ADetails of the QR factorization of the original matrix A as returned by PLASMA_cgeqrf.
[in]TAuxiliary factorization data, computed by PLASMA_cgeqrf.
[in,out]BOn entry, the M-by-N matrix B. On exit, B is overwritten by Q*B or Q**H*B.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_cunmqr
PLASMA_cunmqr_Tile_Async
PLASMA_cunmqr_Tile
PLASMA_dormqr_Tile
PLASMA_sormqr_Tile
PLASMA_cgeqrf_Tile

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.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_cunmqr_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_cunmqr_Tile_Async(side, trans, A, T, B, sequence, &request);
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

Here is the caller graph for this function: