/************************************************************************* * This function performs k-way refinement **************************************************************************/ void Moc_KWayFM(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace, int npasses) { int h, i, ii, iii, j, k, c; int pass, nvtxs, nedges, ncon; int nmoves, nmoved, nswaps, nzgswaps; /* int gnswaps, gnzgswaps; */ int me, firstvtx, lastvtx, yourlastvtx; int from, to = -1, oldto, oldcut, mydomain, yourdomain, imbalanced, overweight; int npes = ctrl->npes, mype = ctrl->mype, nparts = ctrl->nparts; int nlupd, nsupd, nnbrs, nchanged; idxtype *xadj, *ladjncy, *adjwgt, *vtxdist; idxtype *where, *tmp_where, *moved; floattype *lnpwgts, *gnpwgts, *ognpwgts, *pgnpwgts, *movewgts, *overfill; idxtype *update, *supdate, *rupdate, *pe_updates; idxtype *changed, *perm, *pperm, *htable; idxtype *peind, *recvptr, *sendptr; KeyValueType *swchanges, *rwchanges; RInfoType *rinfo, *myrinfo, *tmp_myrinfo, *tmp_rinfo; EdgeType *tmp_edegrees, *my_edegrees, *your_edegrees; floattype lbvec[MAXNCON], *nvwgt, *badmaxpwgt, *ubvec, *tpwgts, lbavg, ubavg; int *nupds_pe; IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->KWayTmr)); /*************************/ /* set up common aliases */ /*************************/ nvtxs = graph->nvtxs; nedges = graph->nedges; ncon = graph->ncon; vtxdist = graph->vtxdist; xadj = graph->xadj; ladjncy = graph->adjncy; adjwgt = graph->adjwgt; firstvtx = vtxdist[mype]; lastvtx = vtxdist[mype+1]; where = graph->where; rinfo = graph->rinfo; lnpwgts = graph->lnpwgts; gnpwgts = graph->gnpwgts; ubvec = ctrl->ubvec; tpwgts = ctrl->tpwgts; nnbrs = graph->nnbrs; peind = graph->peind; recvptr = graph->recvptr; sendptr = graph->sendptr; changed = idxmalloc(nvtxs, "KWR: changed"); rwchanges = wspace->pairs; swchanges = rwchanges + recvptr[nnbrs]; /************************************/ /* set up important data structures */ /************************************/ perm = idxmalloc(nvtxs, "KWR: perm"); pperm = idxmalloc(nparts, "KWR: pperm"); update = idxmalloc(nvtxs, "KWR: update"); supdate = wspace->indices; rupdate = supdate + recvptr[nnbrs]; nupds_pe = imalloc(npes, "KWR: nupds_pe"); htable = idxsmalloc(nvtxs+graph->nrecv, 0, "KWR: lhtable"); badmaxpwgt = fmalloc(nparts*ncon, "badmaxpwgt"); for (i=0; i<nparts; i++) { for (h=0; h<ncon; h++) { badmaxpwgt[i*ncon+h] = ubvec[h]*tpwgts[i*ncon+h]; } } movewgts = fmalloc(nparts*ncon, "KWR: movewgts"); ognpwgts = fmalloc(nparts*ncon, "KWR: ognpwgts"); pgnpwgts = fmalloc(nparts*ncon, "KWR: pgnpwgts"); overfill = fmalloc(nparts*ncon, "KWR: overfill"); moved = idxmalloc(nvtxs, "KWR: moved"); tmp_where = idxmalloc(nvtxs+graph->nrecv, "KWR: tmp_where"); tmp_rinfo = (RInfoType *)GKmalloc(sizeof(RInfoType)*nvtxs, "KWR: tmp_rinfo"); tmp_edegrees = (EdgeType *)GKmalloc(sizeof(EdgeType)*nedges, "KWR: tmp_edegrees"); idxcopy(nvtxs+graph->nrecv, where, tmp_where); for (i=0; i<nvtxs; i++) { tmp_rinfo[i].id = rinfo[i].id; tmp_rinfo[i].ed = rinfo[i].ed; tmp_rinfo[i].ndegrees = rinfo[i].ndegrees; tmp_rinfo[i].degrees = tmp_edegrees+xadj[i]; for (j=0; j<rinfo[i].ndegrees; j++) { tmp_rinfo[i].degrees[j].edge = rinfo[i].degrees[j].edge; tmp_rinfo[i].degrees[j].ewgt = rinfo[i].degrees[j].ewgt; } } nswaps = nzgswaps = 0; /*********************************************************/ /* perform a small number of passes through the vertices */ /*********************************************************/ for (pass=0; pass<npasses; pass++) { if (mype == 0) RandomPermute(nparts, pperm, 1); MPI_Bcast((void *)pperm, nparts, IDX_DATATYPE, 0, ctrl->comm); FastRandomPermute(nvtxs, perm, 1); oldcut = graph->mincut; /* check to see if the partitioning is imbalanced */ Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec); ubavg = savg(ncon, ubvec); lbavg = savg(ncon, lbvec); imbalanced = (lbavg > ubavg) ? 1 : 0; for (c=0; c<2; c++) { scopy(ncon*nparts, gnpwgts, ognpwgts); sset(ncon*nparts, 0.0, movewgts); nmoved = 0; /**********************************************/ /* PASS ONE -- record stats for desired moves */ /**********************************************/ for (iii=0; iii<nvtxs; iii++) { i = perm[iii]; from = tmp_where[i]; nvwgt = graph->nvwgt+i*ncon; for (h=0; h<ncon; h++) if (fabs(nvwgt[h]-gnpwgts[from*ncon+h]) < SMALLFLOAT) break; if (h < ncon) { continue; } /* check for a potential improvement */ if (tmp_rinfo[i].ed >= tmp_rinfo[i].id) { my_edegrees = tmp_rinfo[i].degrees; for (k=0; k<tmp_rinfo[i].ndegrees; k++) { to = my_edegrees[k].edge; if (ProperSide(c, pperm[from], pperm[to])) { for (h=0; h<ncon; h++) if (gnpwgts[to*ncon+h]+nvwgt[h] > badmaxpwgt[to*ncon+h] && nvwgt[h] > 0.0) break; if (h == ncon) break; } } oldto = to; /* check if a subdomain was found that fits */ if (k < tmp_rinfo[i].ndegrees) { for (j=k+1; j<tmp_rinfo[i].ndegrees; j++) { to = my_edegrees[j].edge; if (ProperSide(c, pperm[from], pperm[to])) { for (h=0; h<ncon; h++) if (gnpwgts[to*ncon+h]+nvwgt[h] > badmaxpwgt[to*ncon+h] && nvwgt[h] > 0.0) break; if (h == ncon) { if (my_edegrees[j].ewgt > my_edegrees[k].ewgt || (my_edegrees[j].ewgt == my_edegrees[k].ewgt && IsHBalanceBetterTT(ncon,gnpwgts+oldto*ncon,gnpwgts+to*ncon,nvwgt,ubvec))){ k = j; oldto = my_edegrees[k].edge; } } } } to = oldto; if (my_edegrees[k].ewgt > tmp_rinfo[i].id || (my_edegrees[k].ewgt == tmp_rinfo[i].id && (imbalanced || graph->level > 3 || iii % 8 == 0) && IsHBalanceBetterFT(ncon,gnpwgts+from*ncon,gnpwgts+to*ncon,nvwgt,ubvec))){ /****************************************/ /* Update tmp arrays of the moved vertex */ /****************************************/ tmp_where[i] = to; moved[nmoved++] = i; for (h=0; h<ncon; h++) { lnpwgts[to*ncon+h] += nvwgt[h]; lnpwgts[from*ncon+h] -= nvwgt[h]; gnpwgts[to*ncon+h] += nvwgt[h]; gnpwgts[from*ncon+h] -= nvwgt[h]; movewgts[to*ncon+h] += nvwgt[h]; movewgts[from*ncon+h] -= nvwgt[h]; } tmp_rinfo[i].ed += tmp_rinfo[i].id-my_edegrees[k].ewgt; SWAP(tmp_rinfo[i].id, my_edegrees[k].ewgt, j); if (my_edegrees[k].ewgt == 0) { tmp_rinfo[i].ndegrees--; my_edegrees[k].edge = my_edegrees[tmp_rinfo[i].ndegrees].edge; my_edegrees[k].ewgt = my_edegrees[tmp_rinfo[i].ndegrees].ewgt; } else { my_edegrees[k].edge = from; } /* Update the degrees of adjacent vertices */ for (j=xadj[i]; j<xadj[i+1]; j++) { /* no need to bother about vertices on different pe's */ if (ladjncy[j] >= nvtxs) continue; me = ladjncy[j]; mydomain = tmp_where[me]; myrinfo = tmp_rinfo+me; your_edegrees = myrinfo->degrees; if (mydomain == from) { INC_DEC(myrinfo->ed, myrinfo->id, adjwgt[j]); } else { if (mydomain == to) { INC_DEC(myrinfo->id, myrinfo->ed, adjwgt[j]); } } /* Remove contribution from the .ed of 'from' */ if (mydomain != from) { for (k=0; k<myrinfo->ndegrees; k++) { if (your_edegrees[k].edge == from) { if (your_edegrees[k].ewgt == adjwgt[j]) { myrinfo->ndegrees--; your_edegrees[k].edge = your_edegrees[myrinfo->ndegrees].edge; your_edegrees[k].ewgt = your_edegrees[myrinfo->ndegrees].ewgt; } else { your_edegrees[k].ewgt -= adjwgt[j]; } break; } } } /* Add contribution to the .ed of 'to' */ if (mydomain != to) { for (k=0; k<myrinfo->ndegrees; k++) { if (your_edegrees[k].edge == to) { your_edegrees[k].ewgt += adjwgt[j]; break; } } if (k == myrinfo->ndegrees) { your_edegrees[myrinfo->ndegrees].edge = to; your_edegrees[myrinfo->ndegrees++].ewgt = adjwgt[j]; } } } } } } } /******************************************/ /* Let processors know the subdomain wgts */ /* if all proposed moves commit. */ /******************************************/ MPI_Allreduce((void *)lnpwgts, (void *)pgnpwgts, nparts*ncon, MPI_DOUBLE, MPI_SUM, ctrl->comm); /**************************/ /* compute overfill array */ /**************************/ overweight = 0; for (j=0; j<nparts; j++) { for (h=0; h<ncon; h++) { if (pgnpwgts[j*ncon+h] > ognpwgts[j*ncon+h]) { overfill[j*ncon+h] = (pgnpwgts[j*ncon+h]-badmaxpwgt[j*ncon+h]) / (pgnpwgts[j*ncon+h]-ognpwgts[j*ncon+h]); } else { overfill[j*ncon+h] = 0.0; } overfill[j*ncon+h] = amax(overfill[j*ncon+h], 0.0); overfill[j*ncon+h] *= movewgts[j*ncon+h]; if (overfill[j*ncon+h] > 0.0) overweight = 1; ASSERTP(ctrl, ognpwgts[j*ncon+h] <= badmaxpwgt[j*ncon+h] || pgnpwgts[j*ncon+h] <= ognpwgts[j*ncon+h], (ctrl, "%.4f %.4f %.4f\n", ognpwgts[j*ncon+h], badmaxpwgt[j*ncon+h], pgnpwgts[j*ncon+h])); } } /****************************************************/ /* select moves to undo according to overfill array */ /****************************************************/ if (overweight == 1) { for (iii=0; iii<nmoved; iii++) { i = moved[iii]; oldto = tmp_where[i]; nvwgt = graph->nvwgt+i*ncon; my_edegrees = tmp_rinfo[i].degrees; for (k=0; k<tmp_rinfo[i].ndegrees; k++) if (my_edegrees[k].edge == where[i]) break; for (h=0; h<ncon; h++) if (nvwgt[h] > 0.0 && overfill[oldto*ncon+h] > nvwgt[h]/4.0) break; /**********************************/ /* nullify this move if necessary */ /**********************************/ if (k != tmp_rinfo[i].ndegrees && h != ncon) { moved[iii] = -1; from = oldto; to = where[i]; for (h=0; h<ncon; h++) { overfill[oldto*ncon+h] = amax(overfill[oldto*ncon+h]-nvwgt[h], 0.0); } tmp_where[i] = to; tmp_rinfo[i].ed += tmp_rinfo[i].id-my_edegrees[k].ewgt; SWAP(tmp_rinfo[i].id, my_edegrees[k].ewgt, j); if (my_edegrees[k].ewgt == 0) { tmp_rinfo[i].ndegrees--; my_edegrees[k].edge = my_edegrees[tmp_rinfo[i].ndegrees].edge; my_edegrees[k].ewgt = my_edegrees[tmp_rinfo[i].ndegrees].ewgt; } else { my_edegrees[k].edge = from; } for (h=0; h<ncon; h++) { lnpwgts[to*ncon+h] += nvwgt[h]; lnpwgts[from*ncon+h] -= nvwgt[h]; } /* Update the degrees of adjacent vertices */ for (j=xadj[i]; j<xadj[i+1]; j++) { /* no need to bother about vertices on different pe's */ if (ladjncy[j] >= nvtxs) continue; me = ladjncy[j]; mydomain = tmp_where[me]; myrinfo = tmp_rinfo+me; your_edegrees = myrinfo->degrees; if (mydomain == from) { INC_DEC(myrinfo->ed, myrinfo->id, adjwgt[j]); } else { if (mydomain == to) { INC_DEC(myrinfo->id, myrinfo->ed, adjwgt[j]); } } /* Remove contribution from the .ed of 'from' */ if (mydomain != from) { for (k=0; k<myrinfo->ndegrees; k++) { if (your_edegrees[k].edge == from) { if (your_edegrees[k].ewgt == adjwgt[j]) { myrinfo->ndegrees--; your_edegrees[k].edge = your_edegrees[myrinfo->ndegrees].edge; your_edegrees[k].ewgt = your_edegrees[myrinfo->ndegrees].ewgt; } else { your_edegrees[k].ewgt -= adjwgt[j]; } break; } } } /* Add contribution to the .ed of 'to' */ if (mydomain != to) { for (k=0; k<myrinfo->ndegrees; k++) { if (your_edegrees[k].edge == to) { your_edegrees[k].ewgt += adjwgt[j]; break; } } if (k == myrinfo->ndegrees) { your_edegrees[myrinfo->ndegrees].edge = to; your_edegrees[myrinfo->ndegrees++].ewgt = adjwgt[j]; } } } } } } /*************************************************/ /* PASS TWO -- commit the remainder of the moves */ /*************************************************/ nlupd = nsupd = nmoves = nchanged = 0; for (iii=0; iii<nmoved; iii++) { i = moved[iii]; if (i == -1) continue; where[i] = tmp_where[i]; /* Make sure to update the vertex information */ if (htable[i] == 0) { /* make sure you do the update */ htable[i] = 1; update[nlupd++] = i; } /* Put the vertices adjacent to i into the update array */ for (j=xadj[i]; j<xadj[i+1]; j++) { k = ladjncy[j]; if (htable[k] == 0) { htable[k] = 1; if (k<nvtxs) update[nlupd++] = k; else supdate[nsupd++] = k; } } nmoves++; nswaps++; /* check number of zero-gain moves */ for (k=0; k<rinfo[i].ndegrees; k++) if (rinfo[i].degrees[k].edge == to) break; if (rinfo[i].id == rinfo[i].degrees[k].ewgt) nzgswaps++; if (graph->pexadj[i+1]-graph->pexadj[i] > 0) changed[nchanged++] = i; } /* Tell interested pe's the new where[] info for the interface vertices */ CommChangedInterfaceData(ctrl, graph, nchanged, changed, where, swchanges, rwchanges, wspace->pv4); IFSET(ctrl->dbglvl, DBG_RMOVEINFO, rprintf(ctrl, "\t[%d %d], [%.4f], [%d %d %d]\n", pass, c, badmaxpwgt[0], GlobalSESum(ctrl, nmoves), GlobalSESum(ctrl, nsupd), GlobalSESum(ctrl, nlupd))); /*------------------------------------------------------------- / Time to communicate with processors to send the vertices / whose degrees need to be update. /-------------------------------------------------------------*/ /* Issue the receives first */ for (i=0; i<nnbrs; i++) { MPI_Irecv((void *)(rupdate+sendptr[i]), sendptr[i+1]-sendptr[i], IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->rreq+i); } /* Issue the sends next. This needs some preporcessing */ for (i=0; i<nsupd; i++) { htable[supdate[i]] = 0; supdate[i] = graph->imap[supdate[i]]; } iidxsort(nsupd, supdate); for (j=i=0; i<nnbrs; i++) { yourlastvtx = vtxdist[peind[i]+1]; for (k=j; k<nsupd && supdate[k] < yourlastvtx; k++); MPI_Isend((void *)(supdate+j), k-j, IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->sreq+i); j = k; } /* OK, now get into the loop waiting for the send/recv operations to finish */ MPI_Waitall(nnbrs, ctrl->rreq, ctrl->statuses); for (i=0; i<nnbrs; i++) MPI_Get_count(ctrl->statuses+i, IDX_DATATYPE, nupds_pe+i); MPI_Waitall(nnbrs, ctrl->sreq, ctrl->statuses); /*------------------------------------------------------------- / Place the recieved to-be updated vertices into update[] /-------------------------------------------------------------*/ for (i=0; i<nnbrs; i++) { pe_updates = rupdate+sendptr[i]; for (j=0; j<nupds_pe[i]; j++) { k = pe_updates[j]; if (htable[k-firstvtx] == 0) { htable[k-firstvtx] = 1; update[nlupd++] = k-firstvtx; } } } /*------------------------------------------------------------- / Update the rinfo of the vertices in the update[] array /-------------------------------------------------------------*/ for (ii=0; ii<nlupd; ii++) { i = update[ii]; ASSERT(ctrl, htable[i] == 1); htable[i] = 0; mydomain = where[i]; myrinfo = rinfo+i; tmp_myrinfo = tmp_rinfo+i; my_edegrees = myrinfo->degrees; your_edegrees = tmp_myrinfo->degrees; graph->lmincut -= myrinfo->ed; myrinfo->ndegrees = 0; myrinfo->id = 0; myrinfo->ed = 0; for (j=xadj[i]; j<xadj[i+1]; j++) { yourdomain = where[ladjncy[j]]; if (mydomain != yourdomain) { myrinfo->ed += adjwgt[j]; for (k=0; k<myrinfo->ndegrees; k++) { if (my_edegrees[k].edge == yourdomain) { my_edegrees[k].ewgt += adjwgt[j]; your_edegrees[k].ewgt += adjwgt[j]; break; } } if (k == myrinfo->ndegrees) { my_edegrees[k].edge = yourdomain; my_edegrees[k].ewgt = adjwgt[j]; your_edegrees[k].edge = yourdomain; your_edegrees[k].ewgt = adjwgt[j]; myrinfo->ndegrees++; } ASSERT(ctrl, myrinfo->ndegrees <= xadj[i+1]-xadj[i]); ASSERT(ctrl, tmp_myrinfo->ndegrees <= xadj[i+1]-xadj[i]); } else { myrinfo->id += adjwgt[j]; } } graph->lmincut += myrinfo->ed; tmp_myrinfo->id = myrinfo->id; tmp_myrinfo->ed = myrinfo->ed; tmp_myrinfo->ndegrees = myrinfo->ndegrees; } /* finally, sum-up the partition weights */ MPI_Allreduce((void *)lnpwgts, (void *)gnpwgts, nparts*ncon, MPI_DOUBLE, MPI_SUM, ctrl->comm); } graph->mincut = GlobalSESum(ctrl, graph->lmincut)/2; if (graph->mincut == oldcut) break; } /* gnswaps = GlobalSESum(ctrl, nswaps); gnzgswaps = GlobalSESum(ctrl, nzgswaps); if (mype == 0) printf("niters: %d, nswaps: %d, nzgswaps: %d\n", pass+1, gnswaps, gnzgswaps); */ GKfree((void **)&badmaxpwgt, (void **)&update, (void **)&nupds_pe, (void **)&htable, LTERM); GKfree((void **)&changed, (void **)&pperm, (void **)&perm, (void **)&moved, LTERM); GKfree((void **)&pgnpwgts, (void **)&ognpwgts, (void **)&overfill, (void **)&movewgts, LTERM); GKfree((void **)&tmp_where, (void **)&tmp_rinfo, (void **)&tmp_edegrees, LTERM); IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->KWayTmr)); }
/************************************************************************* * This function performs k-way refinement **************************************************************************/ void MCRandom_KWayEdgeRefineHorizontal(CtrlType *ctrl, GraphType *graph, int nparts, float *orgubvec, int npasses) { int i, ii, iii, j, /*jj,*/ k, /*l,*/ pass, nvtxs, ncon, nmoves, nbnd, myndegrees, same; int from, me, to, oldcut, gain; idxtype *xadj, *adjncy, *adjwgt; idxtype *where, *perm, *bndptr, *bndind; EDegreeType *myedegrees; RInfoType *myrinfo; float *npwgts, *nvwgt, *minwgt, *maxwgt, maxlb, minlb, ubvec[MAXNCON], tvec[MAXNCON]; nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; adjncy = graph->adjncy; adjwgt = graph->adjwgt; bndptr = graph->bndptr; bndind = graph->bndind; where = graph->where; npwgts = graph->npwgts; /* Setup the weight intervals of the various subdomains */ minwgt = fwspacemalloc(ctrl, nparts*ncon); maxwgt = fwspacemalloc(ctrl, nparts*ncon); /* See if the orgubvec consists of identical constraints */ maxlb = minlb = orgubvec[0]; for (i=1; i<ncon; i++) { minlb = (orgubvec[i] < minlb ? orgubvec[i] : minlb); maxlb = (orgubvec[i] > maxlb ? orgubvec[i] : maxlb); } same = (fabs(maxlb-minlb) < .01 ? 1 : 0); /* Let's not get very optimistic. Let Balancing do the work */ ComputeHKWayLoadImbalance(ncon, nparts, npwgts, ubvec); for (i=0; i<ncon; i++) ubvec[i] = amax(ubvec[i], orgubvec[i]); if (!same) { for (i=0; i<nparts; i++) { for (j=0; j<ncon; j++) { maxwgt[i*ncon+j] = ubvec[j]/nparts; minwgt[i*ncon+j] = 1.0/(ubvec[j]*nparts); } } } else { maxlb = ubvec[0]; for (i=1; i<ncon; i++) maxlb = (ubvec[i] > maxlb ? ubvec[i] : maxlb); for (i=0; i<nparts; i++) { for (j=0; j<ncon; j++) { maxwgt[i*ncon+j] = maxlb/nparts; minwgt[i*ncon+j] = 1.0/(maxlb*nparts); } } } perm = idxwspacemalloc(ctrl, nvtxs); if (ctrl->dbglvl&DBG_REFINE) { printf("Partitions: [%5.4f %5.4f], Nv-Nb[%6d %6d]. Cut: %6d, LB: ", npwgts[samin(ncon*nparts, npwgts)], npwgts[samax(ncon*nparts, npwgts)], graph->nvtxs, graph->nbnd, graph->mincut); ComputeHKWayLoadImbalance(ncon, nparts, npwgts, tvec); for (i=0; i<ncon; i++) printf("%.3f ", tvec[i]); printf("\n"); } for (pass=0; pass<npasses; pass++) { ASSERT(ComputeCut(graph, where) == graph->mincut); oldcut = graph->mincut; nbnd = graph->nbnd; RandomPermute(nbnd, perm, 1); for (nmoves=iii=0; iii<graph->nbnd; iii++) { ii = perm[iii]; if (ii >= nbnd) continue; i = bndind[ii]; myrinfo = graph->rinfo+i; if (myrinfo->ed >= myrinfo->id) { /* Total ED is too high */ from = where[i]; nvwgt = graph->nvwgt+i*ncon; if (myrinfo->id > 0 && AreAllHVwgtsBelow(ncon, 1.0, npwgts+from*ncon, -1.0, nvwgt, minwgt+from*ncon)) continue; /* This cannot be moved! */ myedegrees = myrinfo->edegrees; myndegrees = myrinfo->ndegrees; for (k=0; k<myndegrees; k++) { to = myedegrees[k].pid; gain = myedegrees[k].ed - myrinfo->id; if (gain >= 0 && (AreAllHVwgtsBelow(ncon, 1.0, npwgts+to*ncon, 1.0, nvwgt, maxwgt+to*ncon) || IsHBalanceBetterFT(ncon, nparts, npwgts+from*ncon, npwgts+to*ncon, nvwgt, ubvec))) break; } if (k == myndegrees) continue; /* break out if you did not find a candidate */ for (j=k+1; j<myndegrees; j++) { to = myedegrees[j].pid; if ((myedegrees[j].ed > myedegrees[k].ed && (AreAllHVwgtsBelow(ncon, 1.0, npwgts+to*ncon, 1.0, nvwgt, maxwgt+to*ncon) || IsHBalanceBetterFT(ncon, nparts, npwgts+from*ncon, npwgts+to*ncon, nvwgt, ubvec))) || (myedegrees[j].ed == myedegrees[k].ed && IsHBalanceBetterTT(ncon, nparts, npwgts+myedegrees[k].pid*ncon, npwgts+to*ncon, nvwgt, ubvec))) k = j; } to = myedegrees[k].pid; if (myedegrees[k].ed-myrinfo->id == 0 && !IsHBalanceBetterFT(ncon, nparts, npwgts+from*ncon, npwgts+to*ncon, nvwgt, ubvec) && AreAllHVwgtsBelow(ncon, 1.0, npwgts+from*ncon, 0.0, npwgts+from*ncon, maxwgt+from*ncon)) continue; /*===================================================================== * If we got here, we can now move the vertex from 'from' to 'to' *======================================================================*/ graph->mincut -= myedegrees[k].ed-myrinfo->id; IFSET(ctrl->dbglvl, DBG_MOVEINFO, printf("\t\tMoving %6d to %3d. Gain: %4d. Cut: %6d\n", i, to, myedegrees[k].ed-myrinfo->id, graph->mincut)); /* Update where, weight, and ID/ED information of the vertex you moved */ saxpy(ncon, 1.0, nvwgt, 1, npwgts+to*ncon, 1); saxpy(ncon, -1.0, nvwgt, 1, npwgts+from*ncon, 1); where[i] = to; myrinfo->ed += myrinfo->id-myedegrees[k].ed; SWAP(myrinfo->id, myedegrees[k].ed, j); if (myedegrees[k].ed == 0) myedegrees[k] = myedegrees[--myrinfo->ndegrees]; else myedegrees[k].pid = from; if (myrinfo->ed-myrinfo->id < 0) BNDDelete(nbnd, bndind, bndptr, i); /* Update the degrees of adjacent vertices */ for (j=xadj[i]; j<xadj[i+1]; j++) { ii = adjncy[j]; me = where[ii]; myrinfo = graph->rinfo+ii; if (myrinfo->edegrees == NULL) { myrinfo->edegrees = ctrl->wspace.edegrees+ctrl->wspace.cdegree; ctrl->wspace.cdegree += xadj[ii+1]-xadj[ii]; } myedegrees = myrinfo->edegrees; ASSERT(CheckRInfo(myrinfo)); if (me == from) { INC_DEC(myrinfo->ed, myrinfo->id, adjwgt[j]); if (myrinfo->ed-myrinfo->id >= 0 && bndptr[ii] == -1) BNDInsert(nbnd, bndind, bndptr, ii); } else if (me == to) { INC_DEC(myrinfo->id, myrinfo->ed, adjwgt[j]); if (myrinfo->ed-myrinfo->id < 0 && bndptr[ii] != -1) BNDDelete(nbnd, bndind, bndptr, ii); } /* Remove contribution from the .ed of 'from' */ if (me != from) { for (k=0; k<myrinfo->ndegrees; k++) { if (myedegrees[k].pid == from) { if (myedegrees[k].ed == adjwgt[j]) myedegrees[k] = myedegrees[--myrinfo->ndegrees]; else myedegrees[k].ed -= adjwgt[j]; break; } } } /* Add contribution to the .ed of 'to' */ if (me != to) { for (k=0; k<myrinfo->ndegrees; k++) { if (myedegrees[k].pid == to) { myedegrees[k].ed += adjwgt[j]; break; } } if (k == myrinfo->ndegrees) { myedegrees[myrinfo->ndegrees].pid = to; myedegrees[myrinfo->ndegrees++].ed = adjwgt[j]; } } ASSERT(myrinfo->ndegrees <= xadj[ii+1]-xadj[ii]); ASSERT(CheckRInfo(myrinfo)); } nmoves++; } } graph->nbnd = nbnd; if (ctrl->dbglvl&DBG_REFINE) { printf("\t [%5.4f %5.4f], Nb: %6d, Nmoves: %5d, Cut: %6d, LB: ", npwgts[samin(ncon*nparts, npwgts)], npwgts[samax(ncon*nparts, npwgts)], nbnd, nmoves, graph->mincut); ComputeHKWayLoadImbalance(ncon, nparts, npwgts, tvec); for (i=0; i<ncon; i++) printf("%.3f ", tvec[i]); printf("\n"); } if (graph->mincut == oldcut) break; } fwspacefree(ctrl, ncon*nparts); fwspacefree(ctrl, ncon*nparts); idxwspacefree(ctrl, nvtxs); }
/************************************************************************* * This function performs k-way refinement **************************************************************************/ void MCGreedy_KWayEdgeBalanceHorizontal(CtrlType *ctrl, GraphType *graph, int nparts, float *ubvec, int npasses) { int i, ii, /*iii,*/ j, /*jj,*/ k, /*l,*/ pass, nvtxs, ncon, nbnd, myndegrees, oldgain, gain, nmoves; int from, me, to, oldcut; idxtype *xadj, *adjncy, *adjwgt; idxtype *where, *perm, *bndptr, *bndind, *moved; EDegreeType *myedegrees; RInfoType *myrinfo; PQueueType queue; float *npwgts, *nvwgt, *minwgt, *maxwgt, tvec[MAXNCON]; nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; adjncy = graph->adjncy; adjwgt = graph->adjwgt; bndind = graph->bndind; bndptr = graph->bndptr; where = graph->where; npwgts = graph->npwgts; /* Setup the weight intervals of the various subdomains */ minwgt = fwspacemalloc(ctrl, ncon*nparts); maxwgt = fwspacemalloc(ctrl, ncon*nparts); for (i=0; i<nparts; i++) { for (j=0; j<ncon; j++) { maxwgt[i*ncon+j] = ubvec[j]/nparts; minwgt[i*ncon+j] = 1.0/(ubvec[j]*nparts); } } perm = idxwspacemalloc(ctrl, nvtxs); moved = idxwspacemalloc(ctrl, nvtxs); PQueueInit(ctrl, &queue, nvtxs, graph->adjwgtsum[idxamax(nvtxs, graph->adjwgtsum)]); if (ctrl->dbglvl&DBG_REFINE) { printf("Partitions: [%5.4f %5.4f], Nv-Nb[%6d %6d]. Cut: %6d, LB: ", npwgts[samin(ncon*nparts, npwgts)], npwgts[samax(ncon*nparts, npwgts)], graph->nvtxs, graph->nbnd, graph->mincut); ComputeHKWayLoadImbalance(ncon, nparts, npwgts, tvec); for (i=0; i<ncon; i++) printf("%.3f ", tvec[i]); printf("[B]\n"); } for (pass=0; pass<npasses; pass++) { ASSERT(ComputeCut(graph, where) == graph->mincut); /* Check to see if things are out of balance, given the tolerance */ if (MocIsHBalanced(ncon, nparts, npwgts, ubvec)) break; PQueueReset(&queue); idxset(nvtxs, -1, moved); oldcut = graph->mincut; nbnd = graph->nbnd; RandomPermute(nbnd, perm, 1); for (ii=0; ii<nbnd; ii++) { i = bndind[perm[ii]]; PQueueInsert(&queue, i, graph->rinfo[i].ed - graph->rinfo[i].id); moved[i] = 2; } nmoves = 0; for (;;) { if ((i = PQueueGetMax(&queue)) == -1) break; moved[i] = 1; myrinfo = graph->rinfo+i; from = where[i]; nvwgt = graph->nvwgt+i*ncon; if (AreAllHVwgtsBelow(ncon, 1.0, npwgts+from*ncon, -1.0, nvwgt, minwgt+from*ncon)) continue; /* This cannot be moved! */ myedegrees = myrinfo->edegrees; myndegrees = myrinfo->ndegrees; for (k=0; k<myndegrees; k++) { to = myedegrees[k].pid; if (IsHBalanceBetterFT(ncon, nparts, npwgts+from*ncon, npwgts+to*ncon, nvwgt, ubvec)) break; } if (k == myndegrees) continue; /* break out if you did not find a candidate */ for (j=k+1; j<myndegrees; j++) { to = myedegrees[j].pid; if (IsHBalanceBetterTT(ncon, nparts, npwgts+myedegrees[k].pid*ncon, npwgts+to*ncon, nvwgt, ubvec)) k = j; } to = myedegrees[k].pid; j = 0; if (!AreAllHVwgtsBelow(ncon, 1.0, npwgts+from*ncon, 0.0, nvwgt, maxwgt+from*ncon)) j++; if (myedegrees[k].ed-myrinfo->id >= 0) j++; if (!AreAllHVwgtsAbove(ncon, 1.0, npwgts+to*ncon, 0.0, nvwgt, minwgt+to*ncon) && AreAllHVwgtsBelow(ncon, 1.0, npwgts+to*ncon, 1.0, nvwgt, maxwgt+to*ncon)) j++; if (j == 0) continue; /* DELETE if (myedegrees[k].ed-myrinfo->id < 0 && AreAllHVwgtsBelow(ncon, 1.0, npwgts+from*ncon, 0.0, nvwgt, maxwgt+from*ncon) && AreAllHVwgtsAbove(ncon, 1.0, npwgts+to*ncon, 0.0, nvwgt, minwgt+to*ncon) && AreAllHVwgtsBelow(ncon, 1.0, npwgts+to*ncon, 1.0, nvwgt, maxwgt+to*ncon)) continue; */ /*===================================================================== * If we got here, we can now move the vertex from 'from' to 'to' *======================================================================*/ graph->mincut -= myedegrees[k].ed-myrinfo->id; IFSET(ctrl->dbglvl, DBG_MOVEINFO, printf("\t\tMoving %6d to %3d. Gain: %4d. Cut: %6d\n", i, to, myedegrees[k].ed-myrinfo->id, graph->mincut)); /* Update where, weight, and ID/ED information of the vertex you moved */ saxpy(ncon, 1.0, nvwgt, 1, npwgts+to*ncon, 1); saxpy(ncon, -1.0, nvwgt, 1, npwgts+from*ncon, 1); where[i] = to; myrinfo->ed += myrinfo->id-myedegrees[k].ed; SWAP(myrinfo->id, myedegrees[k].ed, j); if (myedegrees[k].ed == 0) myedegrees[k] = myedegrees[--myrinfo->ndegrees]; else myedegrees[k].pid = from; if (myrinfo->ed == 0) BNDDelete(nbnd, bndind, bndptr, i); /* Update the degrees of adjacent vertices */ for (j=xadj[i]; j<xadj[i+1]; j++) { ii = adjncy[j]; me = where[ii]; myrinfo = graph->rinfo+ii; if (myrinfo->edegrees == NULL) { myrinfo->edegrees = ctrl->wspace.edegrees+ctrl->wspace.cdegree; ctrl->wspace.cdegree += xadj[ii+1]-xadj[ii]; } myedegrees = myrinfo->edegrees; ASSERT(CheckRInfo(myrinfo)); oldgain = (myrinfo->ed-myrinfo->id); if (me == from) { INC_DEC(myrinfo->ed, myrinfo->id, adjwgt[j]); if (myrinfo->ed > 0 && bndptr[ii] == -1) BNDInsert(nbnd, bndind, bndptr, ii); } else if (me == to) { INC_DEC(myrinfo->id, myrinfo->ed, adjwgt[j]); if (myrinfo->ed == 0 && bndptr[ii] != -1) BNDDelete(nbnd, bndind, bndptr, ii); } /* Remove contribution from the .ed of 'from' */ if (me != from) { for (k=0; k<myrinfo->ndegrees; k++) { if (myedegrees[k].pid == from) { if (myedegrees[k].ed == adjwgt[j]) myedegrees[k] = myedegrees[--myrinfo->ndegrees]; else myedegrees[k].ed -= adjwgt[j]; break; } } } /* Add contribution to the .ed of 'to' */ if (me != to) { for (k=0; k<myrinfo->ndegrees; k++) { if (myedegrees[k].pid == to) { myedegrees[k].ed += adjwgt[j]; break; } } if (k == myrinfo->ndegrees) { myedegrees[myrinfo->ndegrees].pid = to; myedegrees[myrinfo->ndegrees++].ed = adjwgt[j]; } } /* Update the queue */ if (me == to || me == from) { gain = myrinfo->ed-myrinfo->id; if (moved[ii] == 2) { if (myrinfo->ed > 0) PQueueUpdate(&queue, ii, oldgain, gain); else { PQueueDelete(&queue, ii, oldgain); moved[ii] = -1; } } else if (moved[ii] == -1 && myrinfo->ed > 0) { PQueueInsert(&queue, ii, gain); moved[ii] = 2; } } ASSERT(myrinfo->ndegrees <= xadj[ii+1]-xadj[ii]); ASSERT(CheckRInfo(myrinfo)); } nmoves++; } graph->nbnd = nbnd; if (ctrl->dbglvl&DBG_REFINE) { printf("\t [%5.4f %5.4f], Nb: %6d, Nmoves: %5d, Cut: %6d, LB: ", npwgts[samin(ncon*nparts, npwgts)], npwgts[samax(ncon*nparts, npwgts)], nbnd, nmoves, graph->mincut); ComputeHKWayLoadImbalance(ncon, nparts, npwgts, tvec); for (i=0; i<ncon; i++) printf("%.3f ", tvec[i]); printf("\n"); } if (nmoves == 0) break; } PQueueFree(ctrl, &queue); fwspacefree(ctrl, ncon*nparts); fwspacefree(ctrl, ncon*nparts); idxwspacefree(ctrl, nvtxs); idxwspacefree(ctrl, nvtxs); }
/************************************************************************* * This function performs k-way refinement **************************************************************************/ void Mc_SerialKWayAdaptRefine(GraphType *graph, int nparts, idxtype *home, float *orgubvec, int npasses) { int i, ii, iii, j, k; int nvtxs, ncon, pass, nmoves, myndegrees; int from, me, myhome, to, oldcut, gain, tmp; idxtype *xadj, *adjncy, *adjwgt; idxtype *where; EdgeType *mydegrees; RInfoType *rinfo, *myrinfo; float *npwgts, *nvwgt, *minwgt, *maxwgt, ubvec[MAXNCON]; int gain_is_greater, gain_is_same, fit_in_to, fit_in_from, going_home; int zero_gain, better_balance_ft, better_balance_tt; KeyValueType *cand; int mype; MPI_Comm_rank(MPI_COMM_WORLD, &mype); nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; adjncy = graph->adjncy; adjwgt = graph->adjwgt; where = graph->where; rinfo = graph->rinfo; npwgts = graph->gnpwgts; /* Setup the weight intervals of the various subdomains */ cand = (KeyValueType *)GKmalloc(nvtxs*sizeof(KeyValueType), "cand"); minwgt = fmalloc(nparts*ncon, "minwgt"); maxwgt = fmalloc(nparts*ncon, "maxwgt"); ComputeHKWayLoadImbalance(ncon, nparts, npwgts, ubvec); for (i=0; i<ncon; i++) ubvec[i] = amax(ubvec[i], orgubvec[i]); for (i=0; i<nparts; i++) { for (j=0; j<ncon; j++) { maxwgt[i*ncon+j] = ubvec[j]/(float)nparts; minwgt[i*ncon+j] = ubvec[j]*(float)nparts; } } for (pass=0; pass<npasses; pass++) { oldcut = graph->mincut; for (i=0; i<nvtxs; i++) { cand[i].key = rinfo[i].id-rinfo[i].ed; cand[i].val = i; } ikeysort(nvtxs, cand); nmoves = 0; for (iii=0; iii<nvtxs; iii++) { i = cand[iii].val; myrinfo = rinfo+i; if (myrinfo->ed >= myrinfo->id) { from = where[i]; myhome = home[i]; nvwgt = graph->nvwgt+i*ncon; if (myrinfo->id > 0 && AreAllHVwgtsBelow(ncon, 1.0, npwgts+from*ncon, -1.0, nvwgt, minwgt+from*ncon)) continue; mydegrees = myrinfo->degrees; myndegrees = myrinfo->ndegrees; for (k=0; k<myndegrees; k++) { to = mydegrees[k].edge; gain = mydegrees[k].ewgt - myrinfo->id; if (gain >= 0 && (AreAllHVwgtsBelow(ncon, 1.0, npwgts+to*ncon, 1.0, nvwgt, maxwgt+to*ncon) || IsHBalanceBetterFT(ncon,npwgts+from*ncon,npwgts+to*ncon,nvwgt,ubvec))) { break; } } /* break out if you did not find a candidate */ if (k == myndegrees) continue; for (j=k+1; j<myndegrees; j++) { to = mydegrees[j].edge; going_home = (myhome == to); gain_is_same = (mydegrees[j].ewgt == mydegrees[k].ewgt); gain_is_greater = (mydegrees[j].ewgt > mydegrees[k].ewgt); fit_in_to = AreAllHVwgtsBelow(ncon,1.0,npwgts+to*ncon,1.0,nvwgt,maxwgt+to*ncon); better_balance_ft = IsHBalanceBetterFT(ncon,npwgts+from*ncon, npwgts+to*ncon,nvwgt,ubvec); better_balance_tt = IsHBalanceBetterTT(ncon,npwgts+mydegrees[k].edge*ncon, npwgts+to*ncon,nvwgt,ubvec); if ( (gain_is_greater && (fit_in_to || better_balance_ft) ) || (gain_is_same && ( (fit_in_to && going_home) || better_balance_tt ) ) ) { k = j; } } to = mydegrees[k].edge; going_home = (myhome == to); zero_gain = (mydegrees[k].ewgt == myrinfo->id); fit_in_from = AreAllHVwgtsBelow(ncon,1.0,npwgts+from*ncon,0.0,npwgts+from*ncon, maxwgt+from*ncon); better_balance_ft = IsHBalanceBetterFT(ncon,npwgts+from*ncon, npwgts+to*ncon,nvwgt,ubvec); if (zero_gain && !going_home && !better_balance_ft && fit_in_from) continue; /*===================================================================== * If we got here, we can now move the vertex from 'from' to 'to' *======================================================================*/ graph->mincut -= mydegrees[k].ewgt-myrinfo->id; /* Update where, weight, and ID/ED information of the vertex you moved */ saxpy2(ncon, 1.0, nvwgt, 1, npwgts+to*ncon, 1); saxpy2(ncon, -1.0, nvwgt, 1, npwgts+from*ncon, 1); where[i] = to; myrinfo->ed += myrinfo->id-mydegrees[k].ewgt; SWAP(myrinfo->id, mydegrees[k].ewgt, tmp); if (mydegrees[k].ewgt == 0) { myrinfo->ndegrees--; mydegrees[k].edge = mydegrees[myrinfo->ndegrees].edge; mydegrees[k].ewgt = mydegrees[myrinfo->ndegrees].ewgt; } else mydegrees[k].edge = from; /* Update the degrees of adjacent vertices */ for (j=xadj[i]; j<xadj[i+1]; j++) { ii = adjncy[j]; me = where[ii]; myrinfo = rinfo+ii; mydegrees = myrinfo->degrees; if (me == from) { INC_DEC(myrinfo->ed, myrinfo->id, adjwgt[j]); } else { if (me == to) { INC_DEC(myrinfo->id, myrinfo->ed, adjwgt[j]); } } /* Remove contribution of the ed from 'from' */ if (me != from) { for (k=0; k<myrinfo->ndegrees; k++) { if (mydegrees[k].edge == from) { if (mydegrees[k].ewgt == adjwgt[j]) { myrinfo->ndegrees--; mydegrees[k].edge = mydegrees[myrinfo->ndegrees].edge; mydegrees[k].ewgt = mydegrees[myrinfo->ndegrees].ewgt; } else mydegrees[k].ewgt -= adjwgt[j]; break; } } } /* Add contribution of the ed to 'to' */ if (me != to) { for (k=0; k<myrinfo->ndegrees; k++) { if (mydegrees[k].edge == to) { mydegrees[k].ewgt += adjwgt[j]; break; } } if (k == myrinfo->ndegrees) { mydegrees[myrinfo->ndegrees].edge = to; mydegrees[myrinfo->ndegrees++].ewgt = adjwgt[j]; } } } nmoves++; } } if (graph->mincut == oldcut) break; } GKfree((void **)&minwgt, (void **)&maxwgt, (void **)&cand, LTERM); return; }