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

Functions

void plasma_pcgelqf (plasma_context_t *plasma)
void plasma_pcgelqf_quark (PLASMA_desc A, PLASMA_desc T, 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
Mathieu Faverge
Date:
2010-11-15 c Tue Nov 22 14:35:38 2011

Definition in file pcgelqf.c.


Macro Definition Documentation

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

Definition at line 19 of file pcgelqf.c.

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

Definition at line 20 of file pcgelqf.c.


Function Documentation

void plasma_pcgelqf ( plasma_context_t plasma)

Parallel tile LQ factorization - static scheduling

Definition at line 24 of file pcgelqf.c.

References A, BLKLDD, CORE_cgelqt(), CORE_ctslqt(), CORE_ctsmlq(), CORE_cunmlq(), plasma_desc_t::dtyp, 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_IB, plasma_private_alloc(), plasma_private_free(), PLASMA_RANK, PLASMA_SIZE, PLASMA_SUCCESS, plasma_unpack_args_4, PlasmaConjTrans, PlasmaRight, ss_cond_set, ss_cond_wait, ss_finalize, ss_init, plasma_sequence_t::status, and T.

{
PLASMA_sequence *sequence;
PLASMA_request *request;
int k, m, n;
int next_k;
int next_m;
int next_n;
int ldak, ldam;
int tempkm, tempkn, tempmm, tempnn;
int ib = PLASMA_IB;
PLASMA_Complex32_t *work, *tau;
plasma_unpack_args_4(A, T, sequence, request);
if (sequence->status != PLASMA_SUCCESS)
return;
work = (PLASMA_Complex32_t*)plasma_private_alloc(plasma, ib*T.nb, T.dtyp);
ss_init(A.mt, A.nt, -1);
k = 0;
while (m >= A.mt) {
k++;
m = m-A.mt+k;
}
n = k;
while (k < min(A.mt, A.nt) && m < A.mt) {
next_m = m;
next_n = n;
next_k = k;
next_n++;
if (next_n == A.nt) {
next_m += PLASMA_SIZE;
while (next_m >= A.mt && next_k < min(A.nt, A.mt)) {
next_k++;
next_m = next_m-A.mt+next_k;
}
next_n = next_k;
}
tempkm = k == A.mt-1 ? A.m-k*A.mb : A.mb;
tempkn = k == A.nt-1 ? A.n-k*A.nb : A.nb;
tempmm = m == A.mt-1 ? A.m-m*A.mb : A.mb;
tempnn = n == A.nt-1 ? A.n-n*A.nb : A.nb;
ldak = BLKLDD(A, k);
ldam = BLKLDD(A, m);
if (m == k) {
if (n == k) {
ss_cond_wait(k, k, k-1);
tempkm, tempkn, ib,
A(k, k), ldak,
T(k, k), T.mb,
tau, work);
ss_cond_set(k, k, k);
}
else {
ss_cond_wait(k, n, k-1);
tempkm, tempnn, ib,
A(k, k), ldak,
A(k, n), ldak,
T(k, n), T.mb,
tau, work);
ss_cond_set(k, n, k);
}
}
else {
if (n == k) {
ss_cond_wait(k, k, k);
ss_cond_wait(m, k, k-1);
tempmm, tempkn, tempkn, ib,
A(k, k), ldak,
T(k, k), T.mb,
A(m, k), ldam,
work, T.nb);
}
else {
ss_cond_wait(k, n, k);
ss_cond_wait(m, n, k-1);
tempmm, A.nb, tempmm, tempnn, A.nb, ib,
A(m, k), ldam,
A(m, n), ldam,
A(k, n), ldak,
T(k, n), T.mb,
work, T.nb);
ss_cond_set(m, n, k);
}
}
m = next_m;
n = next_n;
k = next_k;
}
plasma_private_free(plasma, work);
plasma_private_free(plasma, tau);
}

Here is the call graph for this function:

Here is the caller graph for this function:

void plasma_pcgelqf_quark ( PLASMA_desc  A,
PLASMA_desc  T,
PLASMA_sequence sequence,
PLASMA_request request 
)

Parallel tile LQ factorization - dynamic scheduling

Definition at line 137 of file pcgelqf.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, PlasmaRight, plasma_context_struct::quark, QUARK_CORE_cgelqt(), QUARK_CORE_ctslqt(), QUARK_CORE_ctsmlq(), QUARK_CORE_cunmlq(), plasma_sequence_t::quark_sequence, QUARK_Task_Flag_Set(), Quark_Task_Flags_Initializer, plasma_sequence_t::status, T, and TASK_SEQUENCE.

{
int k, m, n;
int ldak, ldam;
int tempkm, tempkn, tempmm, tempnn;
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++) {
tempkm = k == A.mt-1 ? A.m-k*A.mb : A.mb;
tempkn = k == A.nt-1 ? A.n-k*A.nb : A.nb;
ldak = BLKLDD(A, k);
plasma->quark, &task_flags,
tempkm, tempkn, ib, T.nb,
A(k, k), ldak,
T(k, k), T.mb);
for (m = k+1; m < 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, tempkn, ib, T.nb,
A(k, k), ldak,
T(k, k), T.mb,
A(m, k), ldam);
}
for (n = k+1; n < A.nt; n++) {
tempnn = n == A.nt-1 ? A.n-n*A.nb : A.nb;
plasma->quark, &task_flags,
tempkm, tempnn, ib, T.nb,
A(k, k), ldak,
A(k, n), ldak,
T(k, n), T.mb);
for (m = k+1; m < A.mt; m++) {
tempmm = m == A.mt-1 ? A.m-m*A.mb : A.mb;
ldam = BLKLDD(A, m);
plasma->quark, &task_flags,
tempmm, A.nb, tempmm, tempnn, A.mb, ib, T.nb,
A(m, k), ldam,
A(m, n), ldam,
A(k, n), ldak,
T(k, n), T.mb);
}
}
}
}

Here is the call graph for this function:

Here is the caller graph for this function: