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

Go to the source code of this file.

Macros

#define A(m, n)   BLKADDR(A, PLASMA_Complex32_t, (m), (n))
#define T(m, n)   BLKADDR(T, PLASMA_Complex32_t, (m), (n))
#define T2(m, n)   BLKADDR(T, PLASMA_Complex32_t, (m), (n)+A.nt)

Functions

void plasma_pcgeqrfrh_quark (PLASMA_desc A, PLASMA_desc T, int BS, PLASMA_sequence *sequence, PLASMA_request *request)

Detailed Description

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

Version:
2.4.5
Author:
Jakub Kurzak
Hatem Ltaief
Dulceneia Becker
Date:
2010-11-15 c Tue Nov 22 14:35:39 2011

Definition in file pcgeqrfrh.c.


Macro Definition Documentation

#define A (   m,
 
)    BLKADDR(A, PLASMA_Complex32_t, (m), (n))

Definition at line 19 of file pcgeqrfrh.c.

#define T (   m,
 
)    BLKADDR(T, PLASMA_Complex32_t, (m), (n))

Definition at line 20 of file pcgeqrfrh.c.

#define T2 (   m,
 
)    BLKADDR(T, PLASMA_Complex32_t, (m), (n)+A.nt)

Definition at line 21 of file pcgeqrfrh.c.


Function Documentation

void plasma_pcgeqrfrh_quark ( PLASMA_desc  A,
PLASMA_desc  T,
int  BS,
PLASMA_sequence sequence,
PLASMA_request request 
)

Parallel tile QR factorization (reduction Householder) - dynamic scheduling

Definition at line 25 of file pcgeqrfrh.c.

References A, BLKLDD, plasma_desc_t::m, plasma_desc_t::mb, min, plasma_desc_t::mt, plasma_desc_t::n, plasma_desc_t::nb, plasma_desc_t::nt, plasma_context_self(), PLASMA_IB, PLASMA_SUCCESS, PlasmaConjTrans, PlasmaLeft, plasma_context_struct::quark, QUARK_CORE_cgeqrt(), QUARK_CORE_ctsmqr(), QUARK_CORE_ctsqrt(), QUARK_CORE_cttmqr(), QUARK_CORE_cttqrt(), QUARK_CORE_cunmqr(), plasma_sequence_t::quark_sequence, QUARK_Task_Flag_Set(), Quark_Task_Flags_Initializer, plasma_sequence_t::status, T, T2, and TASK_SEQUENCE.

{
int k, m, n;
int M, RD;
int ldaM, ldam, ldaMRD;
int tempkmin, tempkn, tempMm, tempnn, tempmm, tempMRDm;
int ib;
plasma = plasma_context_self();
if (sequence->status != PLASMA_SUCCESS)
return;
QUARK_Task_Flag_Set(&task_flags, TASK_SEQUENCE, (intptr_t)sequence->quark_sequence);
ib = PLASMA_IB;
for (k = 0; k < min(A.mt, A.nt); k++) {
tempkn = k == A.nt-1 ? A.n-k*A.nb : A.nb;
for (M = k; M < A.mt; M += BS) {
tempMm = M == A.mt-1 ? A.m-M*A.mb : A.mb;
tempkmin = min(tempMm, tempkn);
ldaM = BLKLDD(A, M);
plasma->quark, &task_flags,
tempMm, tempkn, ib, T.nb,
A(M, k), ldaM,
T(M, k), T.mb);
for (n = k+1; n < A.nt; n++) {
tempnn = n == A.nt-1 ? A.n-n*A.nb : A.nb;
plasma->quark, &task_flags,
tempMm, tempnn, tempkmin, ib, T.nb,
A(M, k), ldaM,
T(M, k), T.mb,
A(M, n), ldaM);
}
for (m = M+1; m < min(M+BS, A.mt); m++) {
tempmm = m == A.mt-1 ? A.m-m*A.mb : A.mb;
ldam = BLKLDD(A, m);
plasma->quark, &task_flags,
tempmm, tempkn, ib, T.nb,
A(M, k), ldaM,
A(m, k), ldam,
T(m, k), T.mb);
for (n = k+1; n < A.nt; n++) {
tempnn = n == A.nt-1 ? A.n-n*A.nb : A.nb;
plasma->quark, &task_flags,
A.nb, tempnn, tempmm, tempnn, A.nb, ib, T.nb,
A(M, n), ldaM,
A(m, n), ldam,
A(m, k), ldam,
T(m, k), T.mb);
}
}
}
for (RD = BS; RD < A.mt-k; RD *= 2) {
for (M = k; M+RD < A.mt; M += 2*RD) {
tempMRDm = M+RD == A.mt-1 ? A.m-(M+RD)*A.mb : A.mb;
ldaM = BLKLDD(A, M );
ldaMRD = BLKLDD(A, M+RD);
plasma->quark, &task_flags,
tempMRDm, tempkn, ib, T.nb,
A (M , k), ldaM,
A (M+RD, k), ldaMRD,
T2(M+RD, k), T.mb);
for (n = k+1; n < A.nt; n++) {
tempnn = n == A.nt-1 ? A.n-n*A.nb : A.nb;
plasma->quark, &task_flags,
A.nb, tempnn, tempMRDm, tempnn, A.nb, ib, T.nb,
A (M, n), ldaM,
A (M+RD, n), ldaMRD,
A (M+RD, k), ldaMRD,
T2(M+RD, k), T.mb);
}
}
}
}
}

Here is the call graph for this function: