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
core_zgetrf_incpiv.c File Reference
#include <lapacke.h>
#include "common.h"
Include dependency graph for core_zgetrf_incpiv.c:

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)

Detailed Description

PLASMA core_blas kernel PLASMA is a software package provided by Univ. of Tennessee, Univ. of California Berkeley and Univ. of Colorado Denver

Version:
2.4.5
Author:
Hatem Ltaief
Mathieu Faverge
Jakub Kurzak
Date:
2010-11-15 normal z -> c d s

Definition in file core_zgetrf_incpiv.c.


Function Documentation

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.

Parameters:
[in]MThe number of rows of the tile A. M >= 0.
[in]NThe number of columns of the tile A. N >= 0.
[in]IBThe inner-blocking size. IB >= 0.
[in,out]AOn 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]LDAThe leading dimension of the array A. LDA >= max(1,M).
[out]IPIVThe pivot indices; for 1 <= i <= min(M,N), row i of the tile was interchanged with row IPIV(i).
[out]INFOSee returned value.
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if INFO = -k, the k-th argument had an illegal value
>0if 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.

{
int i, j, k, sb;
int iinfo;
/* Check input arguments */
*INFO = 0;
if (M < 0) {
coreblas_error(1, "Illegal value of M");
return -1;
}
if (N < 0) {
coreblas_error(2, "Illegal value of N");
return -2;
}
if (IB < 0) {
coreblas_error(3, "Illegal value of IB");
return -3;
}
if ((LDA < max(1,M)) && (M > 0)) {
coreblas_error(5, "Illegal value of LDA");
return -5;
}
/* Quick return */
if ((M == 0) || (N == 0) || (IB == 0))
k = min(M, N);
for(i =0 ; i < k; i += IB) {
sb = min(IB, k-i);
/*
* Factor diagonal and subdiagonal blocks and test for exact singularity.
*/
iinfo = LAPACKE_zgetf2_work(LAPACK_COL_MAJOR, M-i, sb, &A[LDA*i+i], LDA, &IPIV[i]);
/*
* Adjust INFO and the pivot indices.
*/
if((*INFO == 0) && (iinfo > 0))
*INFO = iinfo + i;
if (i+sb < N) {
M-i, N-(i+sb), sb, sb,
&IPIV[i],
&A[LDA*i+i], LDA,
&A[LDA*(i+sb)+i], LDA);
}
for(j = i; j < i+sb; j++) {
IPIV[j] = i + IPIV[j];
}
}
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

{
int m;
int n;
int ib;
int lda;
int *IPIV;
PLASMA_sequence *sequence;
PLASMA_request *request;
PLASMA_bool check_info;
int iinfo;
int info;
quark_unpack_args_10(quark, m, n, ib, A, lda, IPIV, sequence, request, check_info, iinfo);
CORE_zgetrf_incpiv(m, n, ib, A, lda, IPIV, &info);
if (info != PLASMA_SUCCESS && check_info)
plasma_sequence_flush(quark, sequence, request, iinfo+info);
}

Here is the call graph for this function:

Here is the caller graph for this function:

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.

{
sizeof(int), &m, VALUE,
sizeof(int), &n, VALUE,
sizeof(int), &ib, VALUE,
sizeof(PLASMA_Complex64_t)*nb*nb, A, INOUT,
sizeof(int), &lda, VALUE,
sizeof(int)*nb, IPIV, OUTPUT,
sizeof(PLASMA_sequence*), &sequence, VALUE,
sizeof(PLASMA_request*), &request, VALUE,
sizeof(PLASMA_bool), &check_info, VALUE,
sizeof(int), &iinfo, VALUE,
0);
}

Here is the call graph for this function:

Here is the caller graph for this function: