|
PLASMA
2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
|

Go to the source code of this file.
Functions | |
| int | CORE_zgetrf_incpiv (int M, int N, int IB, PLASMA_Complex64_t *A, int LDA, int *IPIV, int *INFO) |
| void | QUARK_CORE_zgetrf_incpiv (Quark *quark, Quark_Task_Flags *task_flags, int m, int n, int ib, int nb, PLASMA_Complex64_t *A, int lda, int *IPIV, PLASMA_sequence *sequence, PLASMA_request *request, PLASMA_bool check_info, int iinfo) |
| void | CORE_zgetrf_incpiv_quark (Quark *quark) |
PLASMA core_blas kernel PLASMA is a software package provided by Univ. of Tennessee, Univ. of California Berkeley and Univ. of Colorado Denver
Definition in file core_zgetrf_incpiv.c.
| int CORE_zgetrf_incpiv | ( | int | M, |
| int | N, | ||
| int | IB, | ||
| PLASMA_Complex64_t * | A, | ||
| int | LDA, | ||
| int * | IPIV, | ||
| int * | INFO | ||
| ) |
CORE_zgetrf_incpiv computes an LU factorization of a general M-by-N tile A using partial pivoting with row interchanges.
The factorization has the form
A = P * L * U
where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).
This is the right-looking Level 2.5 BLAS version of the algorithm.
| [in] | M | The number of rows of the tile A. M >= 0. |
| [in] | N | The number of columns of the tile A. N >= 0. |
| [in] | IB | The inner-blocking size. IB >= 0. |
| [in,out] | A | On entry, the M-by-N tile to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored. |
| [in] | LDA | The leading dimension of the array A. LDA >= max(1,M). |
| [out] | IPIV | The pivot indices; for 1 <= i <= min(M,N), row i of the tile was interchanged with row IPIV(i). |
| [out] | INFO | See returned value. |
| PLASMA_SUCCESS | successful exit |
| <0 | if INFO = -k, the k-th argument had an illegal value |
| >0 | if INFO = k, U(k,k) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations. |
Definition at line 83 of file core_zgetrf_incpiv.c.
References CORE_zgessm(), coreblas_error, max, min, and PLASMA_SUCCESS.


| void CORE_zgetrf_incpiv_quark | ( | Quark * | quark | ) |
Definition at line 174 of file core_zgetrf_incpiv.c.
References A, CORE_zgetrf_incpiv(), IPIV, plasma_sequence_flush(), PLASMA_SUCCESS, and quark_unpack_args_10.


| void QUARK_CORE_zgetrf_incpiv | ( | Quark * | quark, |
| Quark_Task_Flags * | task_flags, | ||
| int | m, | ||
| int | n, | ||
| int | ib, | ||
| int | nb, | ||
| PLASMA_Complex64_t * | A, | ||
| int | lda, | ||
| int * | IPIV, | ||
| PLASMA_sequence * | sequence, | ||
| PLASMA_request * | request, | ||
| PLASMA_bool | check_info, | ||
| int | iinfo | ||
| ) |
Definition at line 145 of file core_zgetrf_incpiv.c.
References CORE_zgetrf_incpiv_quark(), DAG_CORE_GETRF, INOUT, OUTPUT, QUARK_Insert_Task(), and VALUE.

