5template <
typename Ftype>
21 CHECK_MALLOC(grid3d->iam,
"Enter dAncestorFactor_ASYNC()");
40 donePanelBcast[
i] = 0;
41 donePanelSolve[
i] = 0;
42 localNumChildrenLeft[
i] = 0;
45 for (
int_t k0 = 0; k0 < nnodes; k0++)
47 int_t k = perm_c_supno[k0];
49 int_t ik = myIperm[k_parent];
50 if (ik > -1 && ik < nnodes)
51 localNumChildrenLeft[ik]++;
56 int_t k_st = eTreeTopLims[topoLvl];
57 int_t k_end = eTreeTopLims[topoLvl + 1];
59 for (
int_t k0 = k_st; k0 < k_end; k0++)
61 int_t k = perm_c_supno[k0];
63 dDiagFactorPanelSolve(k, offset, dFBufs);
64 donePanelSolve[k0] = 1;
70 for (
int k0 = k_st; k0 < winSize; ++k0)
72 int_t k = perm_c_supno[k0];
73 int_t offset = k0 % numLA;
74 if (!donePanelBcast[k0])
76 dPanelBcast(k, offset);
77 donePanelBcast[k0] = 1;
83 int_t halfWin = numLA / 2;
88 int_t k = perm_c_supno[k0];
89 int_t offset = (k0 - k1) % winSize;
97 k_upanel = uPanelVec[g2lRow(k)];
99 k_lpanel = lPanelVec[g2lCol(k)];
103 if (UidxSendCounts[k] > 0 && LidxSendCounts[k] > 0)
104 lookAheadUpdate(k, k_parent, k_lpanel, k_upanel);
107 if (k_parent < nsupers)
109 int_t k0_parent = myIperm[k_parent];
110 if (k0_parent > 0 && k0_parent < nnodes)
112 localNumChildrenLeft[k0_parent]--;
113 if (topoLvl < maxTopoLevel - 1 && !localNumChildrenLeft[k0_parent])
116 dDiagFactorPanelSolve(k_parent, dOffset, dFBufs);
117 donePanelSolve[k0_parent] = 1;
123 if (UidxSendCounts[k] > 0 && LidxSendCounts[k] > 0)
124 dSchurCompUpdateExcludeOne(k, k_parent, k_lpanel, k_upanel);
128 for (
int_t k0_next = k1; k0_next <
SUPERLU_MIN(nnodes, k1 + winSize); ++k0_next)
130 int k_next = perm_c_supno[k0_next];
131 if (!localNumChildrenLeft[k0_next])
133 int offset_next = (k0_next - k1) % winSize;
134 if (!(winParity % 2))
135 offset_next += halfWin;
136 dPanelBcast(k_next, offset_next);
137 donePanelBcast[k0_next] = 1;
141 winSize = k0_next - k1;
153 CHECK_MALLOC(grid3d->iam,
"Exit dAncestorFactor_ASYNC()");
Definition: xlupanels.hpp:22
Definition: xlupanels.hpp:176
typename std::conditional< std::is_same< Ftype, double >::value, ddiagFactBufs_t, typename std::conditional< std::is_same< Ftype, float >::value, sdiagFactBufs_t, typename std::conditional< std::is_same< Ftype, doublecomplex >::value, zdiagFactBufs_t, void >::type >::type >::type diagFactBufs_type
Definition: luAuxStructTemplated.hpp:147
Definition: superlu_defs.h:978
int_t * setree
Definition: superlu_defs.h:979
Definition: superlu_defs.h:989
treeTopoInfo_t topoInfo
Definition: superlu_defs.h:999
int_t * nodeList
Definition: superlu_defs.h:992
int_t nNodes
Definition: superlu_defs.h:991
Definition: superlu_defs.h:970
int_t * myIperm
Definition: superlu_defs.h:973
int_t numLvl
Definition: superlu_defs.h:971
int_t * eTreeTopLims
Definition: superlu_defs.h:972
int_t dAncestorFactor(int_t alvl, sForest_t *sforest, diagFactBufs_type< Ftype > **dFBufs, gEtreeInfo_t *gEtreeInfo, int tag_ub)
Definition: dAncestorFactor_impl.hpp:6
Distributed SuperLU data types and function prototypes.
int_t * intMalloc_dist(int_t)
Definition: memory.c:210
int_t getNumLookAhead(superlu_dist_options_t *)
Definition: treeFactorization.c:186
int64_t int_t
Definition: superlu_defs.h:119
int i
Definition: sutil_dist.c:287
#define CHECK_MALLOC(pnum, where)
Definition: util_dist.h:56
#define SUPERLU_MIN(x, y)
Definition: util_dist.h:64
#define SUPERLU_FREE(addr)
Definition: util_dist.h:54