8template <
typename Ftype>
13 int_t myGrid = grid3d->zscp.Iam;
15 int_t sender, receiver;
16 if ((myGrid % (1 << (ilvl + 1))) == 0)
18 sender = myGrid + (1 << ilvl);
24 receiver = myGrid - (1 << ilvl);
28 for (
int_t alvl = ilvl + 1; alvl < maxLvl; ++alvl)
32 int_t numNodes = myNodeCount[alvl];
33 int_t *nodeList = treePerm[alvl];
38 for (
int_t node = 0; node < numNodes; ++node)
40 int_t k0 = nodeList[node];
44 zSendLPanel(k0, receiver);
45 zSendUPanel(k0, receiver);
49 Ftype alpha = one<Ftype>(); Ftype beta = one<Ftype>();
51 zRecvLPanel(k0, sender, alpha, beta);
52 zRecvUPanel(k0, sender, alpha, beta);
61template <
typename Ftype>
65 if (mycol == kcol(k0))
67 int_t lk = g2lCol(k0);
68 if (!lPanelVec[lk].isEmpty())
70 MPI_Send(lPanelVec[lk].blkPtr(0), lPanelVec[lk].nzvalSize(),
71 get_mpi_type<Ftype>(), receiverGrid, k0, grid3d->zscp.comm);
72 SCT->commVolRed += lPanelVec[lk].nzvalSize() *
sizeof(Ftype);
78template <
typename Ftype>
81 if (mycol == kcol(k0))
83 int_t lk = g2lCol(k0);
84 if (!lPanelVec[lk].isEmpty())
88 MPI_Recv(LvalRecvBufs[0], lPanelVec[lk].nzvalSize(), get_mpi_type<Ftype>(), senderGrid, k0,
89 grid3d->zscp.comm, &status);
92 superlu_scal<Ftype>(lPanelVec[lk].nzvalSize(), alpha, lPanelVec[lk].blkPtr(0), 1);
93 superlu_axpy<Ftype>(lPanelVec[lk].nzvalSize(), beta, LvalRecvBufs[0], 1, lPanelVec[lk].blkPtr(0), 1);
99template <
typename Ftype>
103 if (myrow == krow(k0))
105 int_t lk = g2lRow(k0);
106 if (!uPanelVec[lk].isEmpty())
108 MPI_Send(uPanelVec[lk].blkPtr(0), uPanelVec[lk].nzvalSize(),
109 get_mpi_type<Ftype>(), receiverGrid, k0, grid3d->zscp.comm);
110 SCT->commVolRed += uPanelVec[lk].nzvalSize() *
sizeof(Ftype);
116template <
typename Ftype>
119 if (myrow == krow(k0))
121 int_t lk = g2lRow(k0);
122 if (!uPanelVec[lk].isEmpty())
126 MPI_Recv(UvalRecvBufs[0], uPanelVec[lk].nzvalSize(), get_mpi_type<Ftype>(), senderGrid, k0,
127 grid3d->zscp.comm, &status);
130 superlu_scal<Ftype>(uPanelVec[lk].nzvalSize(), alpha, uPanelVec[lk].blkPtr(0), 1);
131 superlu_axpy<Ftype>(uPanelVec[lk].nzvalSize(), beta, UvalRecvBufs[0], 1, uPanelVec[lk].blkPtr(0), 1);
int_t zSendLPanel(int_t k0, int_t receiverGrid)
Definition: lupanels_comm3d_impl.hpp:62
int_t zSendUPanel(int_t k0, int_t receiverGrid)
Definition: lupanels_comm3d_impl.hpp:100
int_t zRecvUPanel(int_t k0, int_t senderGrid, Ftype alpha, Ftype beta)
Definition: lupanels_comm3d_impl.hpp:117
int_t ancestorReduction3d(int_t ilvl, int_t *myNodeCount, int_t **treePerm)
Definition: lupanels_comm3d_impl.hpp:9
int_t zRecvLPanel(int_t k0, int_t senderGrid, Ftype alpha, Ftype beta)
Definition: lupanels_comm3d_impl.hpp:79
Definitions which are precision-neutral.
int_t log2i(int_t index)
Definition: supernodal_etree.c:17
int64_t int_t
Definition: superlu_defs.h:119
double SuperLU_timer_()
Definition: superlu_timer.c:72