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 - Double Real

Functions

int PLASMA_dgebrd_Tile (PLASMA_enum jobu, PLASMA_enum jobvt, PLASMA_desc *A, double *D, double *E, PLASMA_desc *U, PLASMA_desc *VT, PLASMA_desc *T)
int PLASMA_dgelqf_Tile (PLASMA_desc *A, PLASMA_desc *T)
int PLASMA_dgelqs_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_dgels_Tile (PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_dgemm_Tile (PLASMA_enum transA, PLASMA_enum transB, double alpha, PLASMA_desc *A, PLASMA_desc *B, double beta, PLASMA_desc *C)
int PLASMA_dgeqrf_Tile (PLASMA_desc *A, PLASMA_desc *T)
int PLASMA_dgeqrs_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_dgesv_Tile (PLASMA_desc *A, int *IPIV, PLASMA_desc *B)
int PLASMA_dgesv_incpiv_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV, PLASMA_desc *B)
int PLASMA_dgesvd_Tile (PLASMA_enum jobu, PLASMA_enum jobvt, PLASMA_desc *A, double *S, PLASMA_desc *U, PLASMA_desc *VT, PLASMA_desc *T)
int PLASMA_dgetrf_Tile (PLASMA_desc *A, int *IPIV)
int PLASMA_dgetrf_incpiv_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV)
int PLASMA_dgetri_Tile (PLASMA_desc *A, int *IPIV)
int PLASMA_dgetrs_Tile (PLASMA_enum trans, PLASMA_desc *A, int *IPIV, PLASMA_desc *B)
int PLASMA_dgetrs_incpiv_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV, PLASMA_desc *B)
int PLASMA_dlacpy_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
double PLASMA_dlange_Tile (PLASMA_enum norm, PLASMA_desc *A, double *work)
double PLASMA_dlansy_Tile (PLASMA_enum norm, PLASMA_enum uplo, PLASMA_desc *A, double *work)
int PLASMA_dlaset_Tile (PLASMA_enum uplo, double alpha, double beta, PLASMA_desc *A)
int PLASMA_dlaswp_Tile (PLASMA_desc *A, int K1, int K2, int *IPIV, int INCX)
int PLASMA_dlaswpc_Tile (PLASMA_desc *A, int K1, int K2, int *IPIV, int INCX)
int PLASMA_dlauum_Tile (PLASMA_enum uplo, PLASMA_desc *A)
int PLASMA_dorglq_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_dorgqr_Tile (PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_dormlq_Tile (PLASMA_enum side, PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_dormqr_Tile (PLASMA_enum side, PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B)
int PLASMA_dplgsy_Tile (double bump, PLASMA_desc *A, unsigned long long int seed)
int PLASMA_dplrnt_Tile (PLASMA_desc *A, unsigned long long int seed)
int PLASMA_dposv_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_dpotrf_Tile (PLASMA_enum uplo, PLASMA_desc *A)
int PLASMA_dpotri_Tile (PLASMA_enum uplo, PLASMA_desc *A)
int PLASMA_dpotrs_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_dsgesv_Tile (PLASMA_desc *A, int *IPIV, PLASMA_desc *B, PLASMA_desc *X, int *ITER)
int PLASMA_dsposv_Tile (PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B, PLASMA_desc *X, int *ITER)
int PLASMA_dsungesv_Tile (PLASMA_enum trans, PLASMA_desc *A, PLASMA_desc *T, PLASMA_desc *B, PLASMA_desc *X, int *ITER)
int PLASMA_dsyev_Tile (PLASMA_enum jobz, PLASMA_enum uplo, PLASMA_desc *A, double *W, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_dsygst_Tile (PLASMA_enum itype, PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_dsygv_Tile (PLASMA_enum itype, PLASMA_enum jobz, PLASMA_enum uplo, PLASMA_desc *A, PLASMA_desc *B, double *W, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_dsymm_Tile (PLASMA_enum side, PLASMA_enum uplo, double alpha, PLASMA_desc *A, PLASMA_desc *B, double beta, PLASMA_desc *C)
int PLASMA_dsyr2k_Tile (PLASMA_enum uplo, PLASMA_enum trans, double alpha, PLASMA_desc *A, PLASMA_desc *B, double beta, PLASMA_desc *C)
int PLASMA_dsyrk_Tile (PLASMA_enum uplo, PLASMA_enum trans, double alpha, PLASMA_desc *A, double beta, PLASMA_desc *C)
int PLASMA_dsytrd_Tile (PLASMA_enum jobz, PLASMA_enum uplo, PLASMA_desc *A, double *D, double *E, PLASMA_desc *T, PLASMA_desc *Q)
int PLASMA_dtrmm_Tile (PLASMA_enum side, PLASMA_enum uplo, PLASMA_enum transA, PLASMA_enum diag, double alpha, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_dtrsm_Tile (PLASMA_enum side, PLASMA_enum uplo, PLASMA_enum transA, PLASMA_enum diag, double alpha, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_dtrsmpl_Tile (PLASMA_desc *A, PLASMA_desc *L, int *IPIV, PLASMA_desc *B)
int PLASMA_dtrsmrv_Tile (PLASMA_enum side, PLASMA_enum uplo, PLASMA_enum transA, PLASMA_enum diag, double alpha, PLASMA_desc *A, PLASMA_desc *B)
int PLASMA_dtrtri_Tile (PLASMA_enum uplo, PLASMA_enum diag, PLASMA_desc *A)

Detailed Description

This is the group of double real functions using the advanced synchronous interface.


Function Documentation

int PLASMA_dgebrd_Tile ( PLASMA_enum  jobu,
PLASMA_enum  jobvt,
PLASMA_desc A,
double *  D,
double *  E,
PLASMA_desc U,
PLASMA_desc VT,
PLASMA_desc T 
)

PLASMA_dgebrd_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_dgebrd(). 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**T. 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**T (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**T; if JOBVT = 'S', VT contains the first min(m,n) rows of V**T (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_dgebrd
PLASMA_dgebrd_Tile_Async
PLASMA_cgebrd_Tile
PLASMA_dgebrd_Tile
PLASMA_sgebrd_Tile

Definition at line 333 of file dgebrd.c.

References plasma_context_self(), PLASMA_dgebrd_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_dgebrd_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgebrd_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_dgelqf_Tile ( PLASMA_desc A,
PLASMA_desc T 
)

PLASMA_dgelqf_Tile - Computes the tile LQ factorization of a matrix. Tile equivalent of PLASMA_dgelqf(). 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_dgelqs to solve the system of equations.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dgelqf
PLASMA_dgelqf_Tile_Async
PLASMA_cgelqf_Tile
PLASMA_dgelqf_Tile
PLASMA_sgelqf_Tile
PLASMA_dgelqs_Tile

Definition at line 187 of file dgelqf.c.

References plasma_context_self(), PLASMA_dgelqf_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_dgelqf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgelqf_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_dgelqs_Tile ( PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B 
)

PLASMA_dgelqs_Tile - Computes a minimum-norm solution using previously computed LQ factorization. Tile equivalent of PLASMA_dgelqs(). 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_dgelqf.
[in]TAuxiliary factorization data, computed by PLASMA_dgelqf.
[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_dgelqs
PLASMA_dgelqs_Tile_Async
PLASMA_cgelqs_Tile
PLASMA_dgelqs_Tile
PLASMA_sgelqs_Tile
PLASMA_dgelqf_Tile

Definition at line 207 of file dgelqs.c.

References plasma_context_self(), PLASMA_dgelqs_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_dgelqs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgelqs_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_dgels_Tile ( PLASMA_enum  trans,
PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B 
)

PLASMA_dgels_Tile - Solves overdetermined or underdetermined linear system of equations using the tile QR or the tile LQ factorization. Tile equivalent of PLASMA_dgels(). 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; = PlasmaTrans: the linear system involves A**T. 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_dgeqrf; if M < N, A is overwritten by details of its LQ factorization as returned by PLASMA_dgelqf.
[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_dgels
PLASMA_dgels_Tile_Async
PLASMA_cgels_Tile
PLASMA_dgels_Tile
PLASMA_sgels_Tile

Definition at line 267 of file dgels.c.

References plasma_context_self(), PLASMA_dgels_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_dgels_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgels_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_dgemm_Tile ( PLASMA_enum  transA,
PLASMA_enum  transB,
double  alpha,
PLASMA_desc A,
PLASMA_desc B,
double  beta,
PLASMA_desc C 
)

PLASMA_dgemm_Tile - Performs matrix multiplication. Tile equivalent of PLASMA_dgemm(). 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 ugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is transposed; = PlasmaTrans: A is ugate transposed.
[in]transBSpecifies whether the matrix B is transposed, not transposed or ugate transposed: = PlasmaNoTrans: B is not transposed; = PlasmaTrans: B is transposed; = PlasmaTrans: B is ugate 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_dgemm
PLASMA_dgemm_Tile_Async
PLASMA_cgemm_Tile
PLASMA_dgemm_Tile
PLASMA_sgemm_Tile

Definition at line 264 of file dgemm.c.

References plasma_context_self(), PLASMA_dgemm_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_dgemm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgemm_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_dgeqrf_Tile ( PLASMA_desc A,
PLASMA_desc T 
)

PLASMA_dgeqrf_Tile - Computes the tile QR factorization of a matrix. Tile equivalent of PLASMA_dgeqrf(). 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_dgeqrs to solve the system of equations.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dgeqrf
PLASMA_dgeqrf_Tile_Async
PLASMA_cgeqrf_Tile
PLASMA_dgeqrf_Tile
PLASMA_sgeqrf_Tile
PLASMA_dgeqrs_Tile

Definition at line 186 of file dgeqrf.c.

References plasma_context_self(), PLASMA_dgeqrf_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_dgeqrf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgeqrf_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_dgeqrs_Tile ( PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B 
)

PLASMA_dgeqrs_Tile - Computes a minimum-norm solution using the tile QR factorization. Tile equivalent of PLASMA_dgeqrf(). 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_dgeqrf.
[in]TAuxiliary factorization data, computed by PLASMA_dgeqrf.
[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_dgeqrs
PLASMA_dgeqrs_Tile_Async
PLASMA_cgeqrs_Tile
PLASMA_dgeqrs_Tile
PLASMA_sgeqrs_Tile
PLASMA_dgeqrf_Tile

Definition at line 206 of file dgeqrs.c.

References plasma_context_self(), PLASMA_dgeqrs_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_dgeqrs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgeqrs_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_dgesv_incpiv_Tile ( PLASMA_desc A,
PLASMA_desc L,
int *  IPIV,
PLASMA_desc B 
)

PLASMA_dgesv_incpiv_Tile - Solves a system of linear equations using the tile LU factorization. Tile equivalent of PLASMA_dgetrf_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_dgesv_incpiv
PLASMA_dgesv_incpiv_Tile_Async
PLASMA_cgesv_incpiv_Tile
PLASMA_dgesv_incpiv_Tile
PLASMA_sgesv_incpiv_Tile
PLASMA_dcgesv_Tile

Definition at line 203 of file dgesv_incpiv.c.

References plasma_context_self(), PLASMA_dgesv_incpiv_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_dgesv_incpiv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgesv_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_dgesv_Tile ( PLASMA_desc A,
int *  IPIV,
PLASMA_desc B 
)

PLASMA_dgesv_Tile - Solves a system of linear equations using the tile LU factorization. Tile equivalent of PLASMA_dgetrf(). 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_dgesv
PLASMA_dgesv_Tile_Async
PLASMA_cgesv_Tile
PLASMA_dgesv_Tile
PLASMA_sgesv_Tile
PLASMA_dcgesv_Tile

Definition at line 187 of file dgesv.c.

References plasma_context_self(), PLASMA_dgesv_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_dgesv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgesv_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_dgesvd_Tile ( PLASMA_enum  jobu,
PLASMA_enum  jobvt,
PLASMA_desc A,
double *  S,
PLASMA_desc U,
PLASMA_desc VT,
PLASMA_desc T 
)

PLASMA_dgesvd_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_dgesvd(). 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**T. 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**T (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**T; if JOBVT = 'S', VT contains the first min(m,n) rows of V**T (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_dgesvd
PLASMA_dgesvd_Tile_Async
PLASMA_cgesvd_Tile
PLASMA_dgesvd_Tile
PLASMA_sgesvd_Tile

Definition at line 333 of file dgesvd.c.

References plasma_context_self(), PLASMA_dgesvd_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_dgesvd_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgesvd_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_dgetrf_incpiv_Tile ( PLASMA_desc A,
PLASMA_desc L,
int *  IPIV 
)

PLASMA_dgetrf_incpiv_Tile - Computes the tile LU factorization of a matrix. Tile equivalent of PLASMA_dgetrf_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_dgetrs_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_dgetrf_incpiv
PLASMA_dgetrf_incpiv_Tile_Async
PLASMA_cgetrf_incpiv_Tile
PLASMA_dgetrf_incpiv_Tile
PLASMA_sgetrf_incpiv_Tile
PLASMA_dgetrs_incpiv_Tile

Definition at line 184 of file dgetrf_incpiv.c.

References plasma_context_self(), PLASMA_dgetrf_incpiv_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_dgetrf_incpiv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgetrf_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_dgetrf_Tile ( PLASMA_desc A,
int *  IPIV 
)

PLASMA_dgetrf_Tile - Computes the tile LU factorization of a matrix. Tile equivalent of PLASMA_dgetrf(). 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_dgetrs 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_dgetrf
PLASMA_dgetrf_Tile_Async
PLASMA_cgetrf_Tile
PLASMA_dgetrf_Tile
PLASMA_sgetrf_Tile
PLASMA_dgetrs_Tile

Definition at line 189 of file dgetrf.c.

References plasma_context_self(), PLASMA_dgetrf_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_dgetrf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgetrf_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_dgetri_Tile ( PLASMA_desc A,
int *  IPIV 
)

PLASMA_dgetri_Tile - Computes the inverse of a matrix using the LU factorization computed by PLASMA_dgetrf. 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_dgetri(). 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_dgetrf. 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_dgetrf.
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_dgetri
PLASMA_dgetri_Tile_Async
PLASMA_cgetri_Tile
PLASMA_dgetri_Tile
PLASMA_sgetri_Tile
PLASMA_dgetrf_Tile

Definition at line 175 of file dgetri.c.

References PLASMA_Alloc_Workspace_dgetri_Tile_Async(), plasma_context_self(), plasma_desc_mat_free(), PLASMA_dgetri_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;
PLASMA_desc descW;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_dgetri_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
/* Allocate workspace */
PLASMA_dgetri_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_dgetrs_incpiv_Tile ( PLASMA_desc A,
PLASMA_desc L,
int *  IPIV,
PLASMA_desc B 
)

PLASMA_dgetrs_incpiv_Tile - Solves a system of linear equations using previously computed LU factorization. Tile equivalent of PLASMA_dgetrs_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_dgetrf_incpiv.
[in]LAuxiliary factorization data, related to the tile L factor, computed by PLASMA_dgetrf_incpiv.
[in]IPIVThe pivot indices from PLASMA_dgetrf_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_dgetrs_incpiv
PLASMA_dgetrs_incpiv_Tile_Async
PLASMA_cgetrs_incpiv_Tile
PLASMA_dgetrs_incpiv_Tile
PLASMA_sgetrs_incpiv_Tile
PLASMA_dgetrf_incpiv_Tile

Definition at line 206 of file dgetrs_incpiv.c.

References plasma_context_self(), PLASMA_dgetrs_incpiv_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_dgetrs_incpiv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgetrs_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_dgetrs_Tile ( PLASMA_enum  trans,
PLASMA_desc A,
int *  IPIV,
PLASMA_desc B 
)

PLASMA_dgetrs_Tile - Solves a system of linear equations using previously computed LU factorization. Tile equivalent of PLASMA_dgetrs(). 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) = PlasmaTrans: A**T * X = B (Conjugate transpose)
[in]AThe tile factors L and U from the factorization, computed by PLASMA_dgetrf.
[in]IPIVThe pivot indices from PLASMA_dgetrf.
[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_dgetrs
PLASMA_dgetrs_Tile_Async
PLASMA_cgetrs_Tile
PLASMA_dgetrs_Tile
PLASMA_sgetrs_Tile
PLASMA_dgetrf_Tile

Definition at line 199 of file dgetrs.c.

References plasma_context_self(), PLASMA_dgetrs_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_dgetrs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dgetrs_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_dlacpy_Tile ( PLASMA_enum  uplo,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dlacpy_Tile - Tile equivalent of PLASMA_dlacpy(). 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_dlacpy
PLASMA_dlacpy_Tile_Async
PLASMA_clacpy_Tile
PLASMA_dlacpy_Tile
PLASMA_slacpy_Tile

Definition at line 183 of file dlacpy.c.

References plasma_context_self(), PLASMA_dlacpy_Tile_Async(), 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_dlacpy_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dlacpy_Tile_Async(uplo, A, B, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
}

Here is the call graph for this function:

double PLASMA_dlange_Tile ( PLASMA_enum  norm,
PLASMA_desc A,
double *  work 
)

PLASMA_dlange_Tile - Tile equivalent of PLASMA_dlange(). 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]workdouble precision 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_dlange
PLASMA_dlange_Tile_Async
PLASMA_clange_Tile
PLASMA_dlange_Tile
PLASMA_slange_Tile

Definition at line 193 of file dlange.c.

References plasma_context_self(), PLASMA_dlange_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), and plasma_sequence_destroy().

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

Here is the call graph for this function:

double PLASMA_dlansy_Tile ( PLASMA_enum  norm,
PLASMA_enum  uplo,
PLASMA_desc A,
double *  work 
)

PLASMA_dlansy_Tile - Tile equivalent of PLASMA_dlansy(). 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]workdouble precision array of dimension PLASMA_SIZE is PLASMA_STATIC_SCHEDULING is used, and NULL otherwise.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dlansy
PLASMA_dlansy_Tile_Async
PLASMA_clansy_Tile
PLASMA_dlansy_Tile
PLASMA_slansy_Tile

Definition at line 195 of file dlansy.c.

References plasma_context_self(), PLASMA_dlansy_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), and plasma_sequence_destroy().

{
PLASMA_sequence *sequence = NULL;
double value;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_dlansy_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dlansy_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_dlaset_Tile ( PLASMA_enum  uplo,
double  alpha,
double  beta,
PLASMA_desc A 
)

PLASMA_dlaset_Tile - Tile equivalent of PLASMA_dlaset(). 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_dlaset
PLASMA_dlaset_Tile_Async
PLASMA_claset_Tile
PLASMA_dlaset_Tile
PLASMA_slaset_Tile

Definition at line 172 of file dlaset.c.

References plasma_context_self(), PLASMA_dlaset_Tile_Async(), 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_dlaset_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dlaset_Tile_Async(uplo, alpha, beta, A, sequence, &request);
plasma_sequence_destroy(plasma, sequence);
}

Here is the call graph for this function:

int PLASMA_dlaswp_Tile ( PLASMA_desc A,
int  K1,
int  K2,
int *  IPIV,
int  INCX 
)

PLASMA_dlaswp_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_dlaswp(). 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_dgetrf.
[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_dgetrf.
[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_dlaswp
PLASMA_dlaswp_Tile_Async
PLASMA_claswp_Tile
PLASMA_dlaswp_Tile
PLASMA_slaswp_Tile
PLASMA_dgetrf_Tile

Definition at line 176 of file dlaswp.c.

References plasma_context_self(), PLASMA_dlaswp_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_dlaswp_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dlaswp_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_dlaswpc_Tile ( PLASMA_desc A,
int  K1,
int  K2,
int *  IPIV,
int  INCX 
)

PLASMA_dlaswpc_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_dlaswpc(). 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_dgetrf.
[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_dgetrf.
[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_dlaswpc
PLASMA_dlaswpc_Tile_Async
PLASMA_claswpc_Tile
PLASMA_dlaswpc_Tile
PLASMA_slaswpc_Tile
PLASMA_dgetrf_Tile

Definition at line 176 of file dlaswpc.c.

References plasma_context_self(), PLASMA_dlaswpc_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_dlaswpc_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dlaswpc_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_dlauum_Tile ( PLASMA_enum  uplo,
PLASMA_desc A 
)

PLASMA_dlauum_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_dlauum(). 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_dlauum
PLASMA_dlauum_Tile_Async
PLASMA_clauum_Tile
PLASMA_dlauum_Tile
PLASMA_slauum_Tile
PLASMA_dpotri_Tile

Definition at line 172 of file dlauum.c.

References plasma_context_self(), PLASMA_dlauum_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_dlauum_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dlauum_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_dorglq_Tile ( PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B 
)

PLASMA_dorglq_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_dgelqf. 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_dgelqf.
[in]TAuxiliary factorization data, computed by PLASMA_dgelqf.
[out]BOn exit, the M-by-N matrix Q.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dorglq
PLASMA_dorglq_Tile_Async
PLASMA_cunglq_Tile
PLASMA_dorglq_Tile
PLASMA_sorglq_Tile
PLASMA_dgelqf_Tile

Definition at line 202 of file dorglq.c.

References plasma_context_self(), PLASMA_dorglq_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_dorglq_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dorglq_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_dorgqr_Tile ( PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc Q 
)

PLASMA_dorgqr_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_dgeqrf. 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_dgeqrf.
[in]TAuxiliary factorization data, computed by PLASMA_dgeqrf.
[out]QOn exit, the M-by-N matrix Q.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dorgqr
PLASMA_dorgqr_Tile_Async
PLASMA_cungqr_Tile
PLASMA_dorgqr_Tile
PLASMA_sorgqr_Tile
PLASMA_dgeqrf_Tile

Definition at line 200 of file dorgqr.c.

References plasma_context_self(), PLASMA_dorgqr_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_dorgqr_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dorgqr_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_dormlq_Tile ( PLASMA_enum  side,
PLASMA_enum  trans,
PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B 
)

PLASMA_dormlq_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_dgelqf_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**T from the left; = PlasmaRight: apply Q or Q**T from the right. Currently only PlasmaLeft is supported.
[in]transIntended usage: = PlasmaNoTrans: no transpose, apply Q; = PlasmaTrans: ugate transpose, apply Q**T. Currently only PlasmaTrans is supported.
[in]ADetails of the LQ factorization of the original matrix A as returned by PLASMA_dgelqf.
[in]TAuxiliary factorization data, computed by PLASMA_dgelqf.
[in,out]BOn entry, the M-by-N matrix B. On exit, B is overwritten by Q*B or Q**T*B.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dormlq
PLASMA_dormlq_Tile_Async
PLASMA_cunmlq_Tile
PLASMA_dormlq_Tile
PLASMA_sormlq_Tile
PLASMA_dgelqf_Tile

Definition at line 247 of file dormlq.c.

References plasma_context_self(), PLASMA_dormlq_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_dormlq_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dormlq_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_dormqr_Tile ( PLASMA_enum  side,
PLASMA_enum  trans,
PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B 
)

PLASMA_dormqr_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_dgeqrf_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**T from the left; = PlasmaRight: apply Q or Q**T from the right. Currently only PlasmaLeft is supported.
[in]transIntended usage: = PlasmaNoTrans: no transpose, apply Q; = PlasmaTrans: ugate transpose, apply Q**T. Currently only PlasmaTrans is supported.
[in]ADetails of the QR factorization of the original matrix A as returned by PLASMA_dgeqrf.
[in]TAuxiliary factorization data, computed by PLASMA_dgeqrf.
[in,out]BOn entry, the M-by-N matrix B. On exit, B is overwritten by Q*B or Q**T*B.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dormqr
PLASMA_dormqr_Tile_Async
PLASMA_cormqr_Tile
PLASMA_dormqr_Tile
PLASMA_sormqr_Tile
PLASMA_dgeqrf_Tile

Definition at line 250 of file dormqr.c.

References plasma_context_self(), PLASMA_dormqr_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_dormqr_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dormqr_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_dplgsy_Tile ( double  bump,
PLASMA_desc A,
unsigned long long int  seed 
)

PLASMA_dplgsy_Tile - Generate a random hermitian matrix by tiles. Tile equivalent of PLASMA_dplgsy(). 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_dplgsy
PLASMA_dplgsy_Tile_Async
PLASMA_cplgsy_Tile
PLASMA_dplgsy_Tile
PLASMA_splgsy_Tile
PLASMA_dplrnt_Tile
PLASMA_dplgsy_Tile

Definition at line 153 of file dplgsy.c.

References plasma_context_self(), PLASMA_dplgsy_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_dplgsy_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dplgsy_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_dplrnt_Tile ( PLASMA_desc A,
unsigned long long int  seed 
)

PLASMA_dplrnt_Tile - Generate a random matrix by tiles. Tile equivalent of PLASMA_dplrnt(). 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_dplrnt
PLASMA_dplrnt_Tile_Async
PLASMA_cplrnt_Tile
PLASMA_dplrnt_Tile
PLASMA_splrnt_Tile
PLASMA_dplgsy_Tile
PLASMA_dplgsy_Tile

Definition at line 151 of file dplrnt.c.

References plasma_context_self(), PLASMA_dplrnt_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_dplrnt_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dplrnt_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_dposv_Tile ( PLASMA_enum  uplo,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dposv_Tile - Solves a symmetric positive definite or Hermitian positive definite system of linear equations using the Cholesky factorization. Tile equivalent of PLASMA_dposv(). 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**T*U or A = L*L**T.
[in,out]BOn entry, the N-by-NRHS right hand side matrix B. On exit, if return value = 0, the N-by-NRHS solution matrix X.
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_dposv
PLASMA_dposv_Tile_Async
PLASMA_cposv_Tile
PLASMA_dposv_Tile
PLASMA_sposv_Tile

Definition at line 212 of file dposv.c.

References plasma_context_self(), PLASMA_dposv_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_dposv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dposv_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_dpotrf_Tile ( PLASMA_enum  uplo,
PLASMA_desc A 
)

PLASMA_dpotrf_Tile - Computes the Cholesky factorization of a symmetric positive definite or Hermitian positive definite matrix. Tile equivalent of PLASMA_dpotrf(). 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**T*U or A = L*L**T.
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_dpotrf
PLASMA_dpotrf_Tile_Async
PLASMA_cpotrf_Tile
PLASMA_dpotrf_Tile
PLASMA_spotrf_Tile
PLASMA_dpotrs_Tile

Definition at line 183 of file dpotrf.c.

References plasma_context_self(), PLASMA_dpotrf_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_dpotrf_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dpotrf_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_dpotri_Tile ( PLASMA_enum  uplo,
PLASMA_desc A 
)

PLASMA_dpotri_Tile - Computes the inverse of a complex Hermitian positive definite matrix A using the Cholesky factorization A = U**T*U or A = L*L**T computed by PLASMA_dpotrf. Tile equivalent of PLASMA_dpotri(). 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**T*U or A = L*L**T, as computed by PLASMA_dpotrf. 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_dpotri
PLASMA_dpotri_Tile_Async
PLASMA_cpotri_Tile
PLASMA_dpotri_Tile
PLASMA_spotri_Tile
PLASMA_dpotrf_Tile

Definition at line 172 of file dpotri.c.

References plasma_context_self(), PLASMA_dpotri_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_dpotri_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dpotri_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_dpotrs_Tile ( PLASMA_enum  uplo,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dpotrs_Tile - Solves a system of linear equations using previously computed Cholesky factorization. Tile equivalent of PLASMA_dpotrs(). 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**T*U or A = L*L**T, computed by PLASMA_dpotrf.
[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_dpotrs
PLASMA_dpotrs_Tile_Async
PLASMA_cpotrs_Tile
PLASMA_dpotrs_Tile
PLASMA_spotrs_Tile
PLASMA_dpotrf_Tile

Definition at line 187 of file dpotrs.c.

References plasma_context_self(), PLASMA_dpotrs_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_dpotrs_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dpotrs_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_dsgesv_Tile ( PLASMA_desc A,
int *  IPIV,
PLASMA_desc B,
PLASMA_desc X,
int *  ITER 
)

PLASMA_dsgesv_Tile - Solves a system of linear equations using the tile LU factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_dsgesv(). 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.
  • If the iterative refinement converged, A is not modified;
  • otherwise, it fell back to double precision solution, and then A contains the tile L and U factors from the factorization (not equivalent to LAPACK).
[out]LOn exit:
  • if the iterative refinement converged, L is not modified;
  • otherwise, it fell back to double precision solution, and then L is an auxiliary factorization data, related to the tile L factor, necessary to solve the system of equations (not equivalent to LAPACK).
[out]IPIVOn exit, the pivot indices that define the permutations (not equivalent to LAPACK).
[in]BOn entry, the N-by-NRHS matrix of right hand side matrix B.
[out]XOn exit, if return value = 0, the N-by-NRHS solution matrix X.
[out]ITERThe number of the current iteration in the iterative refinement process
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_dsgesv
PLASMA_dsgesv_Tile_Async
PLASMA_dsgesv_Tile
PLASMA_dgesv_Tile

Definition at line 374 of file dsgesv.c.

References plasma_context_self(), PLASMA_dsgesv_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, and plasma_sequence_t::status.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_dsgesv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
status = PLASMA_dsgesv_Tile_Async(A, IPIV, B, X, ITER, sequence, &request);
if (status != PLASMA_SUCCESS)
return status;
status = sequence->status;
plasma_sequence_destroy(plasma, sequence);
return status;
}

Here is the call graph for this function:

int PLASMA_dsposv_Tile ( PLASMA_enum  uplo,
PLASMA_desc A,
PLASMA_desc B,
PLASMA_desc X,
int *  ITER 
)

PLASMA_dsposv_Tile - Solves a symmetric positive definite or Hermitian positive definite system of linear equations using the Cholesky factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_dsposv(). 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 N-by-N symmetric positive definite (or Hermitian) coefficient matrix A. If uplo = PlasmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced.
  • If the iterative refinement converged, A is not modified;
  • otherwise, it falled backed to double precision solution,
[in]BOn entry, the N-by-NRHS matrix of right hand side matrix B.
[out]XOn exit, if return value = 0, the N-by-NRHS solution matrix X.
[out]ITERThe number of the current iteration in the iterative refinement process
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_dsposv
PLASMA_dsposv_Tile_Async
PLASMA_dsposv_Tile
PLASMA_dposv_Tile

Definition at line 323 of file dsposv.c.

References plasma_context_self(), PLASMA_dsposv_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, and plasma_sequence_t::status.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_dsposv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
status = PLASMA_dsposv_Tile_Async(uplo, A, B, X, ITER, sequence, &request);
if (status != PLASMA_SUCCESS)
return status;
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_dsungesv_Tile ( PLASMA_enum  trans,
PLASMA_desc A,
PLASMA_desc T,
PLASMA_desc B,
PLASMA_desc X,
int *  ITER 
)

PLASMA_dsungesv_Tile - Solves symmetric linear system of equations using the tile QR or the tile LQ factorization and mixed-precision iterative refinement. Tile equivalent of PLASMA_dsungesv(). 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; = PlasmaTrans: the linear system involves A**H. Currently only PlasmaNoTrans is supported.
[in,out]A
  • If the iterative refinement converged, A is not modified;
  • otherwise, it fell back to double precision solution, and on exit the M-by-N matrix A contains: if M >= N, A is overwritten by details of its QR factorization as returned by PLASMA_dgeqrf; if M < N, A is overwritten by details of its LQ factorization as returned by PLASMA_dgelqf.
[out]TOn exit:
  • if the iterative refinement converged, T is not modified;
  • otherwise, it fell back to double precision solution, and then T is an auxiliary factorization data.
[in,out]BOn entry, the M-by-NRHS matrix B of right hand side vectors, stored columnwise;
[in]BThe N-by-NRHS matrix of right hand side matrix B.
[out]XIf return value = 0, X is the solution vectors, stored columnwise: if M >= N, rows 1 to N of X contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of the modulus of elements N+1 to M in that column; if M < N, rows 1 to N of X contain the minimum norm solution vectors;
[out]ITERThe number of the current iteration in the iterative refinement process
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
See also:
PLASMA_dsungesv
PLASMA_dsungesv_Tile_Async
PLASMA_dsungesv_Tile
PLASMA_dgels_Tile

Definition at line 328 of file dsungesv.c.

References plasma_context_self(), PLASMA_dsungesv_Tile_Async(), plasma_dynamic_sync, PLASMA_ERR_NOT_INITIALIZED, plasma_fatal_error(), PLASMA_REQUEST_INITIALIZER, plasma_sequence_create(), plasma_sequence_destroy(), PLASMA_SUCCESS, and plasma_sequence_t::status.

{
PLASMA_sequence *sequence = NULL;
int status;
plasma = plasma_context_self();
if (plasma == NULL) {
plasma_fatal_error("PLASMA_dsungesv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
status = PLASMA_dsungesv_Tile_Async(trans, A, T, B, X, ITER, sequence, &request);
if (status != PLASMA_SUCCESS)
return status;
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_dsyev_Tile ( PLASMA_enum  jobz,
PLASMA_enum  uplo,
PLASMA_desc A,
double *  W,
PLASMA_desc T,
PLASMA_desc Q 
)

PLASMA_dsyev_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_dsyev 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_dsyev_Tile
PLASMA_dsyev_Tile_Async
PLASMA_cheev
PLASMA_dsyev
PLASMA_ssyev

Definition at line 274 of file dsyev.c.

References plasma_context_self(), PLASMA_dsyev_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_dsyev_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsyev_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_dsygst_Tile ( PLASMA_enum  itype,
PLASMA_enum  uplo,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dsygst_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**T)*A*inv(U) or inv(L)*A*inv(L**T) 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**T or L**T*A*L. B must have been previously factorized as U**T*U or L*L**T by PLASMA_DPOTRF. ONLY PlasmaItype == 1 and PlasmaLower supported! Tile equivalent of PLASMA_dsygst(). 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_DPOTRF.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value
See also:
PLASMA_dsygst
PLASMA_dsygst_Tile_Async
PLASMA_chegst_Tile
PLASMA_dsygst_Tile
PLASMA_ssygst_Tile
PLASMA_dsygst_Tile

Definition at line 233 of file dsygst.c.

References plasma_context_self(), PLASMA_dsygst_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_dsygst_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsygst_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_dsygv_Tile ( PLASMA_enum  itype,
PLASMA_enum  jobz,
PLASMA_enum  uplo,
PLASMA_desc A,
PLASMA_desc B,
double *  W,
PLASMA_desc T,
PLASMA_desc Q 
)

PLASMA_dsygv_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_dsygv(). 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**T*B*Z = I; if ITYPE = 3, Z**T*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**T*U or B = L*L**T.
[in]LDBThe leading dimension of the array B. LDA >= max(1,N).
[in,out]TOn entry, descriptor as return by PLASMA_Alloc_Workspace_dsygv 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_dsygv 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_dsygv
PLASMA_dsygv_Tile_Async
PLASMA_chegv_Tile
PLASMA_dsygv_Tile
PLASMA_ssygv_Tile

Definition at line 369 of file dsygv.c.

References plasma_context_self(), PLASMA_dsygv_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_dsygv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsygv_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_dsymm_Tile ( PLASMA_enum  side,
PLASMA_enum  uplo,
double  alpha,
PLASMA_desc A,
PLASMA_desc B,
double  beta,
PLASMA_desc C 
)

PLASMA_dsymm_Tile - Performs symmetric matrix multiplication. Tile equivalent of PLASMA_dsymm(). 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_dsymm
PLASMA_dsymm_Tile_Async
PLASMA_csymm_Tile
PLASMA_dsymm_Tile
PLASMA_ssymm_Tile

Definition at line 254 of file dsymm.c.

References plasma_context_self(), PLASMA_dsymm_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_dsymm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsymm_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_dsyr2k_Tile ( PLASMA_enum  uplo,
PLASMA_enum  trans,
double  alpha,
PLASMA_desc A,
PLASMA_desc B,
double  beta,
PLASMA_desc C 
)

PLASMA_dsyr2k_Tile - Performs symmetric rank k update. Tile equivalent of PLASMA_dsyr2k(). 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 ugate transposed: = PlasmaNoTrans: A is not transposed; = PlasmaTrans: A is ugate 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_dsyr2k_Tile
PLASMA_csyr2k
PLASMA_dsyr2k
PLASMA_ssyr2k

Definition at line 250 of file dsyr2k.c.

References plasma_context_self(), PLASMA_dsyr2k_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_dsyr2k_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsyr2k_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_dsyrk_Tile ( PLASMA_enum  uplo,
PLASMA_enum  trans,
double  alpha,
PLASMA_desc A,
double  beta,
PLASMA_desc C 
)

PLASMA_dsyrk_Tile - Performs rank k update. Tile equivalent of PLASMA_dsyrk(). 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 ugate 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_dsyrk_Tile
PLASMA_csyrk
PLASMA_dsyrk
PLASMA_ssyrk

Definition at line 227 of file dsyrk.c.

References plasma_context_self(), PLASMA_dsyrk_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_dsyrk_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsyrk_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_dsytrd_Tile ( PLASMA_enum  jobz,
PLASMA_enum  uplo,
PLASMA_desc A,
double *  D,
double *  E,
PLASMA_desc T,
PLASMA_desc Q 
)

PLASMA_dsytrd_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**T * A * Q = S. Tile equivalent of PLASMA_dsytrd(). 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_dsytrd
PLASMA_dsytrd_Tile_Async
PLASMA_chetrd_Tile
PLASMA_dsytrd_Tile
PLASMA_ssytrd_Tile
PLASMA_dsyev_Tile

Definition at line 279 of file dsytrd.c.

References plasma_context_self(), PLASMA_dsytrd_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_dsytrd_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dsytrd_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_dtrmm_Tile ( PLASMA_enum  side,
PLASMA_enum  uplo,
PLASMA_enum  transA,
PLASMA_enum  diag,
double  alpha,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dtrmm_Tile - Computes triangular solve. Tile equivalent of PLASMA_dtrmm(). 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 ugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaTrans: A is ugate 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_dtrmm
PLASMA_dtrmm_Tile_Async
PLASMA_ctrmm_Tile
PLASMA_dtrmm_Tile
PLASMA_strmm_Tile

Definition at line 249 of file dtrmm.c.

References plasma_context_self(), PLASMA_dtrmm_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_dtrmm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dtrmm_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_dtrsm_Tile ( PLASMA_enum  side,
PLASMA_enum  uplo,
PLASMA_enum  transA,
PLASMA_enum  diag,
double  alpha,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dtrsm_Tile - Computes triangular solve. Tile equivalent of PLASMA_dtrsm(). 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 ugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaTrans: A is ugate 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_dtrsm
PLASMA_dtrsm_Tile_Async
PLASMA_ctrsm_Tile
PLASMA_dtrsm_Tile
PLASMA_strsm_Tile

Definition at line 249 of file dtrsm.c.

References plasma_context_self(), PLASMA_dtrsm_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_dtrsm_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dtrsm_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_dtrsmpl_Tile ( PLASMA_desc A,
PLASMA_desc L,
int *  IPIV,
PLASMA_desc B 
)

PLASMA_dtrsmpl_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_dgetrf_incpiv.
[in]LAuxiliary factorization data, related to the tile L factor, computed by PLASMA_dgetrf_incpiv.
[in]IPIVThe pivot indices from PLASMA_dgetrf_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_dtrsmpl
PLASMA_dtrsmpl_Tile_Async
PLASMA_ctrsmpl_Tile
PLASMA_dtrsmpl_Tile
PLASMA_strsmpl_Tile
PLASMA_dgetrf_incpiv_Tile

Definition at line 191 of file dtrsmpl.c.

References plasma_context_self(), PLASMA_dtrsmpl_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_dtrsmpl_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dtrsmpl_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_dtrsmrv_Tile ( PLASMA_enum  side,
PLASMA_enum  uplo,
PLASMA_enum  transA,
PLASMA_enum  diag,
double  alpha,
PLASMA_desc A,
PLASMA_desc B 
)

PLASMA_dtrsmrv_Tile - Computes triangular solve. Tile equivalent of PLASMA_dtrsmrv(). 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 ugate transposed: = PlasmaNoTrans: A is transposed; = PlasmaTrans: A is not transposed; = PlasmaTrans: A is ugate 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_dtrsmrv
PLASMA_dtrsmrv_Tile_Async
PLASMA_ctrsmrv_Tile
PLASMA_dtrsmrv_Tile
PLASMA_strsmrv_Tile

Definition at line 249 of file dtrsmrv.c.

References plasma_context_self(), PLASMA_dtrsmrv_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_dtrsmrv_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dtrsmrv_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_dtrtri_Tile ( PLASMA_enum  uplo,
PLASMA_enum  diag,
PLASMA_desc A 
)

PLASMA_dtrtri_Tile - Computes the inverse of a complex upper or lower triangular matrix A. Tile equivalent of PLASMA_dtrtri(). 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_dtrtri
PLASMA_dtrtri_Tile_Async
PLASMA_ctrtri_Tile
PLASMA_dtrtri_Tile
PLASMA_strtri_Tile
PLASMA_dpotri_Tile

Definition at line 191 of file dtrtri.c.

References plasma_context_self(), PLASMA_dtrtri_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_dtrtri_Tile", "PLASMA not initialized");
}
plasma_sequence_create(plasma, &sequence);
PLASMA_dtrtri_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: