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

Go to the source code of this file.
Functions | |
| int | CORE_dtstrf (int M, int N, int IB, int NB, double *U, int LDU, double *A, int LDA, double *L, int LDL, int *IPIV, double *WORK, int LDWORK, int *INFO) |
| void | QUARK_CORE_dtstrf (Quark *quark, Quark_Task_Flags *task_flags, int m, int n, int ib, int nb, double *U, int ldu, double *A, int lda, double *L, int ldl, int *IPIV, PLASMA_sequence *sequence, PLASMA_request *request, PLASMA_bool check_info, int iinfo) |
| void | CORE_dtstrf_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_dtstrf.c.
| int CORE_dtstrf | ( | int | M, |
| int | N, | ||
| int | IB, | ||
| int | NB, | ||
| double * | U, | ||
| int | LDU, | ||
| double * | A, | ||
| int | LDA, | ||
| double * | L, | ||
| int | LDL, | ||
| int * | IPIV, | ||
| double * | WORK, | ||
| int | LDWORK, | ||
| int * | INFO | ||
| ) |
CORE_dtstrf computes an LU factorization of a complex matrix formed by an upper triangular NB-by-N tile U on top of a M-by-N tile A using partial pivoting with row interchanges.
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] | NB | |
| [in,out] | U | On entry, the NB-by-N upper triangular tile. On exit, the new factor U from the factorization |
| [in] | LDU | The leading dimension of the array U. LDU >= max(1,NB). |
| [in,out] | A | On entry, the M-by-N tile to be factored. On exit, the factor L from the factorization |
| [in] | LDA | The leading dimension of the array A. LDA >= max(1,M). |
| [in,out] | L | On entry, the IB-by-N lower triangular tile. On exit, the interchanged rows form the tile A in case of pivoting. |
| [in] | LDL | The leading dimension of the array L. LDL >= max(1,IB). |
| [out] | IPIV | The pivot indices; for 1 <= i <= min(M,N), row i of the tile U was interchanged with row IPIV(i) of the tile A. |
| [in,out] | WORK | |
| [in] | LDWORK | The dimension of the array WORK. |
| [out] | INFO |
| 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 99 of file core_dtstrf.c.
References cblas_dcopy(), cblas_dger(), cblas_dscal(), cblas_dswap(), cblas_idamax(), CblasColMajor, CORE_dssssm(), coreblas_error, max, min, and PLASMA_SUCCESS.


| void CORE_dtstrf_quark | ( | Quark * | quark | ) |
Definition at line 258 of file core_dtstrf.c.
References A, CORE_dtstrf(), IPIV, L, plasma_sequence_flush(), PLASMA_SUCCESS, and quark_unpack_args_17.


| void QUARK_CORE_dtstrf | ( | Quark * | quark, |
| Quark_Task_Flags * | task_flags, | ||
| int | m, | ||
| int | n, | ||
| int | ib, | ||
| int | nb, | ||
| double * | U, | ||
| int | ldu, | ||
| double * | A, | ||
| int | lda, | ||
| double * | L, | ||
| int | ldl, | ||
| int * | IPIV, | ||
| PLASMA_sequence * | sequence, | ||
| PLASMA_request * | request, | ||
| PLASMA_bool | check_info, | ||
| int | iinfo | ||
| ) |
Definition at line 220 of file core_dtstrf.c.
References CORE_dtstrf_quark(), DAG_CORE_TSTRF, INOUT, LOCALITY, OUTPUT, QUARK_Insert_Task(), QUARK_REGION_D, QUARK_REGION_U, SCRATCH, and VALUE.

