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

Go to the source code of this file.

Macros

#define COMPLEX

Functions

int CORE_zttlqt (int M, int N, int IB, PLASMA_Complex64_t *A1, int LDA1, PLASMA_Complex64_t *A2, int LDA2, PLASMA_Complex64_t *T, int LDT, PLASMA_Complex64_t *TAU, PLASMA_Complex64_t *WORK)
void QUARK_CORE_zttlqt (Quark *quark, Quark_Task_Flags *task_flags, int m, int n, int ib, int nb, PLASMA_Complex64_t *A1, int lda1, PLASMA_Complex64_t *A2, int lda2, PLASMA_Complex64_t *T, int ldt)
void CORE_zttlqt_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
Dulceneia Becker
Date:
2010-11-15 normal z -> c d s

Definition in file core_zttlqt.c.


Macro Definition Documentation

#define COMPLEX

Definition at line 20 of file core_zttlqt.c.


Function Documentation

int CORE_zttlqt ( int  M,
int  N,
int  IB,
PLASMA_Complex64_t A1,
int  LDA1,
PLASMA_Complex64_t A2,
int  LDA2,
PLASMA_Complex64_t T,
int  LDT,
PLASMA_Complex64_t TAU,
PLASMA_Complex64_t WORK 
)

CORE_zttlqt computes a LQ factorization of a rectangular matrix formed by coupling side-by-side a complex M-by-M lower triangular tile A1 and a complex M-by-N lower triangular tile A2:

| A1 A2 | = L * Q

The tile Q is represented as a product of elementary reflectors

Q = H(k)' . . . H(2)' H(1)', 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; conjg(v(i+1:n)) is stored on exit in A2(i,1:n), and tau in TAU(i).

Parameters:
[in]MThe number of rows of the tile A1 and A2. M >= 0. The number of columns of the tile A1.
[in]NThe number of columns of the tile A2. N >= 0.
[in]IBThe inner-blocking size. IB >= 0.
[in,out]A1On entry, the M-by-M tile A1. On exit, the elements on and below the diagonal of the array contain the M-by-M lower trapezoidal tile L; the elements above the diagonal are not referenced.
[in]LDA1The leading dimension of the array A1. LDA1 >= max(1,N).
[in,out]A2On entry, the M-by-N lower triangular tile A2. On exit, the elements on and below the diagonal of the array with the array TAU, represent the unitary tile Q as a product of elementary reflectors (see Further Details).
[in]LDA2The leading dimension of the array A2. LDA2 >= 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).
[in,out]WORK
Returns:
Return values:
PLASMA_SUCCESSsuccessful exit
<0if -i, the i-th argument had an illegal value

Definition at line 100 of file core_zttlqt.c.

References CBLAS_SADDR, cblas_zaxpy(), cblas_zcopy(), cblas_zgemv(), cblas_zgerc(), cblas_ztrmv(), CblasColMajor, CORE_zlaset(), CORE_zparfb(), CORE_zpemv(), coreblas_error, max, min, PLASMA_SUCCESS, PlasmaForward, PlasmaNonUnit, PlasmaNoTrans, PlasmaRight, PlasmaRowwise, PlasmaUpper, and PlasmaUpperLower.

{
static PLASMA_Complex64_t zone = 1.0;
static PLASMA_Complex64_t zzero = 0.0;
#ifdef COMPLEX
static int ione = 1;
#endif
int i, j, l, ii, sb, mi, ni;
/* 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) {
coreblas_error(3, "Illegal value of IB");
return -3;
}
if ((LDA2 < max(1,M)) && (M > 0)) {
coreblas_error(7, "Illegal value of LDA2");
return -7;
}
/* Quick return */
if ((M == 0) || (N == 0) || (IB == 0))
/* TODO: Need to check why some cases require
* this to not have uninitialized values */
0., 0., T, LDT);
for(ii = 0; ii < M; ii += IB) {
sb = min(M-ii, IB);
for(i = 0; i < sb; i++) {
j = ii + i;
mi = sb-i-1;
ni = min( j + 1, N);
/*
* Generate elementary reflector H( II*IB+I ) to annihilate A( II*IB+I, II*IB+I:M ).
*/
#ifdef COMPLEX
LAPACKE_zlacgv_work(ni, &A2[j], LDA2);
LAPACKE_zlacgv_work(ione, &A1[LDA1*j+j], LDA1);
#endif
LAPACKE_zlarfg_work(ni+1, &A1[LDA1*j+j], &A2[j], LDA2, &TAU[j]);
if (mi > 0) {
/*
* Apply H( j-1 ) to A( j:II+IB-1, j-1:M ) from the right.
*/
mi,
&A1[LDA1*j+(j+1)], 1,
WORK, 1);
mi, ni,
CBLAS_SADDR(zone), &A2[j+1], LDA2,
&A2[j], LDA2,
CBLAS_SADDR(zone), WORK, 1);
alpha = -(TAU[j]);
mi, CBLAS_SADDR(alpha),
WORK, 1,
&A1[LDA1*j+j+1], 1);
CblasColMajor, mi, ni,
CBLAS_SADDR(alpha), WORK, 1,
&A2[j], LDA2,
&A2[j+1], LDA2);
}
/*
* Calculate T.
*/
if (i > 0 ) {
l = min(i, max(0, N-ii));
alpha = -(TAU[j]);
PlasmaNoTrans, PlasmaRowwise,
i , min(j, N), l,
alpha, &A2[ii], LDA2,
&A2[j], LDA2,
zzero, &T[LDT*j], 1,
WORK);
/* T(0:i-1, j) = T(0:i-1, ii:j-1) * T(0:i-1, j) */
(CBLAS_TRANSPOSE)PlasmaNoTrans,
i, &T[LDT*ii], LDT,
&T[LDT*j], 1);
}
#ifdef COMPLEX
LAPACKE_zlacgv_work(ni, &A2[j], LDA2 );
LAPACKE_zlacgv_work(ione, &A1[LDA1*j+j], LDA1 );
#endif
T[LDT*j+i] = TAU[j];
}
/* Apply Q to the rest of the matrix to the right */
if (M > ii+sb) {
mi = M-(ii+sb);
ni = min(ii+sb, N);
l = min(sb, max(0, ni-ii));
PlasmaRight, PlasmaNoTrans,
mi, IB, mi, ni, sb, l,
&A1[LDA1*ii+ii+sb], LDA1,
&A2[ii+sb], LDA2,
&A2[ii], LDA2,
&T[LDT*ii], LDT,
WORK, M);
}
}
}

Here is the call graph for this function:

Here is the caller graph for this function:

void CORE_zttlqt_quark ( Quark quark)

Definition at line 273 of file core_zttlqt.c.

References CORE_zttlqt(), quark_unpack_args_11, T, and TAU.

{
int m;
int n;
int ib;
int lda1;
int lda2;
int ldt;
quark_unpack_args_11(quark, m, n, ib, A1, lda1, A2, lda2, T, ldt, TAU, WORK);
CORE_zttlqt(m, n, ib, A1, lda1, A2, lda2, T, ldt, TAU, WORK);
}

Here is the call graph for this function:

Here is the caller graph for this function:

void QUARK_CORE_zttlqt ( Quark quark,
Quark_Task_Flags task_flags,
int  m,
int  n,
int  ib,
int  nb,
PLASMA_Complex64_t A1,
int  lda1,
PLASMA_Complex64_t A2,
int  lda2,
PLASMA_Complex64_t T,
int  ldt 
)

Definition at line 244 of file core_zttlqt.c.

References CORE_zttlqt_quark(), DAG_CORE_TTLQT, INOUT, LOCALITY, OUTPUT, QUARK_Insert_Task(), QUARK_REGION_D, QUARK_REGION_L, SCRATCH, and VALUE.

{
sizeof(int), &m, VALUE,
sizeof(int), &n, VALUE,
sizeof(int), &ib, VALUE,
sizeof(int), &lda1, VALUE,
sizeof(int), &lda2, VALUE,
sizeof(PLASMA_Complex64_t)*ib*nb, T, OUTPUT,
sizeof(int), &ldt, VALUE,
sizeof(PLASMA_Complex64_t)*nb, NULL, SCRATCH,
sizeof(PLASMA_Complex64_t)*ib*nb, NULL, SCRATCH,
0);
}

Here is the call graph for this function:

Here is the caller graph for this function: