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

Go to the source code of this file.

Functions

int CORE_sgeqrt (int M, int N, int IB, float *A, int LDA, float *T, int LDT, float *TAU, float *WORK)
void QUARK_CORE_sgeqrt (Quark *quark, Quark_Task_Flags *task_flags, int m, int n, int ib, int nb, float *A, int lda, float *T, int ldt)
void CORE_sgeqrt_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 s Tue Nov 22 14:35:16 2011

Definition in file core_sgeqrt.c.


Function Documentation

int CORE_sgeqrt ( int  M,
int  N,
int  IB,
float *  A,
int  LDA,
float *  T,
int  LDT,
float *  TAU,
float *  WORK 
)

CORE_sgeqrt computes a QR factorization of a complex M-by-N tile A: A = Q * R.

The tile Q is represented as a product of elementary reflectors

Q = H(1) H(2) . . . H(k), where k = min(M,N).

Each H(i) has the form

H(i) = I - tau * v * v'

where tau is a complex scalar, and v is a complex vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).

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 A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal tile R (R is upper triangular if M >= N); the elements below the diagonal, with the array TAU, represent the unitary tile Q as a product of elementary reflectors (see Further Details).
[in]LDAThe leading dimension of the array A. LDA >= max(1,M).
[out]TThe IB-by-N triangular factor T of the block reflector. T is upper triangular by block (economic storage); The rest of the array is not referenced.
[in]LDTThe leading dimension of the array T. LDT >= IB.
[out]TAUThe scalar factors of the elementary reflectors (see Further Details).
[out]WORK
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value

Definition at line 86 of file core_sgeqrt.c.

References coreblas_error, lapack_const, max, min, PLASMA_SUCCESS, PlasmaColumnwise, PlasmaForward, PlasmaLeft, and PlasmaTrans.

{
int i, k, sb;
/* Check input arguments */
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) || ( (IB == 0) && ((M > 0) && (N > 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;
}
if ((LDT < max(1,IB)) && (IB > 0)) {
coreblas_error(7, "Illegal value of LDT");
return -7;
}
/* 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);
LAPACKE_sgeqr2_work(LAPACK_COL_MAJOR, M-i, sb,
&A[LDA*i+i], LDA, &TAU[i], WORK);
LAPACKE_slarft_work(LAPACK_COL_MAJOR,
M-i, sb,
&A[LDA*i+i], LDA, &TAU[i],
&T[LDT*i], LDT);
if (N > i+sb) {
LAPACKE_slarfb_work(
LAPACK_COL_MAJOR,
M-i, N-i-sb, sb,
&A[LDA*i+i], LDA,
&T[LDT*i], LDT,
&A[LDA*(i+sb)+i], LDA,
WORK, N-i-sb);
}
}
}

Here is the caller graph for this function:

void CORE_sgeqrt_quark ( Quark quark)

Definition at line 181 of file core_sgeqrt.c.

References A, CORE_sgeqrt(), quark_unpack_args_9, T, and TAU.

{
int m;
int n;
int ib;
float *A;
int lda;
float *T;
int ldt;
float *TAU;
float *WORK;
quark_unpack_args_9(quark, m, n, ib, A, lda, T, ldt, TAU, WORK);
CORE_sgeqrt(m, n, ib, A, lda, T, ldt, TAU, WORK);
}

Here is the call graph for this function:

Here is the caller graph for this function:

void QUARK_CORE_sgeqrt ( Quark quark,
Quark_Task_Flags task_flags,
int  m,
int  n,
int  ib,
int  nb,
float *  A,
int  lda,
float *  T,
int  ldt 
)

Definition at line 155 of file core_sgeqrt.c.

References CORE_sgeqrt_quark(), DAG_CORE_GEQRT, INOUT, OUTPUT, QUARK_Insert_Task(), SCRATCH, and VALUE.

{
sizeof(int), &m, VALUE,
sizeof(int), &n, VALUE,
sizeof(int), &ib, VALUE,
sizeof(float)*nb*nb, A, INOUT,
sizeof(int), &lda, VALUE,
sizeof(float)*ib*nb, T, OUTPUT,
sizeof(int), &ldt, VALUE,
sizeof(float)*nb, NULL, SCRATCH,
sizeof(float)*ib*nb, NULL, SCRATCH,
0);
}

Here is the call graph for this function:

Here is the caller graph for this function: