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

Go to the source code of this file.

Macros

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

Functions

void plasma_pcungqr_quark (PLASMA_desc A, PLASMA_desc Q, 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:
Hatem Ltaief
Jakub Kurzak
Mathieu Faverge
Date:
2010-11-15 c Tue Nov 22 14:35:41 2011

Definition in file pcungqr.c.


Macro Definition Documentation

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

Definition at line 19 of file pcungqr.c.

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

Definition at line 20 of file pcungqr.c.

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

Definition at line 21 of file pcungqr.c.


Function Documentation

void plasma_pcungqr_quark ( PLASMA_desc  A,
PLASMA_desc  Q,
PLASMA_desc  T,
PLASMA_sequence sequence,
PLASMA_request request 
)

Parallel construction of Q using tile V (application to identity) - dynamic scheduling

Definition at line 25 of file pcungqr.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, PlasmaLeft, PlasmaNoTrans, Q, plasma_context_struct::quark, QUARK_CORE_ctsmqr(), QUARK_CORE_cunmqr(), 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, ldqk, ldam, ldqm;
int tempmm, tempnn, tempkmin, tempkm;
int tempAkm, tempAkn;
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 = min(A.mt, A.nt)-1; k >= 0; k--) {
tempAkm = k == A.mt-1 ? A.m-k*A.mb : A.mb;
tempAkn = k == A.nt-1 ? A.n-k*A.nb : A.nb;
tempkmin = min( tempAkn, tempAkm );
tempkm = k == Q.mt-1 ? Q.m-k*Q.mb : Q.mb;
ldak = BLKLDD(A, k);
ldqk = BLKLDD(Q, k);
for (m = Q.mt - 1; m > k; m--) {
tempmm = m == Q.mt-1 ? Q.m-m*Q.mb : Q.mb;
ldam = BLKLDD(A, m);
ldqm = BLKLDD(Q, m);
for (n = 0; n < Q.nt; n++) {
tempnn = n == Q.nt-1 ? Q.n-n*Q.nb : Q.nb;
plasma->quark, &task_flags,
Q.mb, tempnn, tempmm, tempnn, tempAkn, ib, T.nb,
Q(k, n), ldqk,
Q(m, n), ldqm,
A(m, k), ldam,
T(m, k), T.mb);
}
}
for (n = 0; n < Q.nt; n++) {
tempnn = n == Q.nt-1 ? Q.n-n*Q.nb : Q.nb;
plasma->quark, &task_flags,
tempkm, tempnn, tempkmin, ib, T.nb,
A(k, k), ldak,
T(k, k), T.mb,
Q(k, n), ldqk);
}
}
}

Here is the call graph for this function: