/*********************************************************************************** * This function determines the cost of moving data between the two meshes assuming * that a good matching between the two partitions was done! ************************************************************************************/ int ComputeMapCost(idxtype nvtxs, idxtype nparts, idxtype *fepart, idxtype *cpart) { idxtype i, j, k, n, ncomm; KeyValueType cand[nparts*nparts]; idxtype fmatched[nparts], cmatched[nparts]; /* Compute the overlap */ for (i=0; i<nparts; i++) { for (j=0; j<nparts; j++) { cand[i*nparts+j].key = 0; cand[i*nparts+j].val = i*nparts+j; } } for (k=0, i=0; i<nvtxs; i++) { if (cpart[i] >= 0) { cand[(fepart[i]-1)*nparts+(cpart[i]-1)].key++; k++; } } mprintf("Contact points: %D\n", k); ikeysort(nparts*nparts, cand); iset(nparts, -1, fmatched); iset(nparts, -1, cmatched); for (ncomm=0, k=nparts*nparts-1; k>=0; k--) { i = cand[k].val/nparts; j = cand[k].val%nparts; if (fmatched[i] == -1 && cmatched[j] == -1) { fmatched[i] = j; cmatched[j] = i; } else ncomm += cand[k].key; } mprintf("Ncomm: %D\n", ncomm); return ncomm; }
/************************************************************************* * This function computes the subdomain graph **************************************************************************/ void EliminateSubDomainEdges(CtrlType *ctrl, GraphType *graph, int nparts, float *tpwgts) { int i, ii, j, k, me, other, nvtxs, total, max, avg, totalout, nind, ncand, ncand2, target, target2, nadd; int min, move, cpwgt, tvwgt; idxtype *xadj, *adjncy, *vwgt, *adjwgt, *pwgts, *where, *maxpwgt, *pmat, *ndoms, *mypmat, *otherpmat, *ind; KeyValueType *cand, *cand2; nvtxs = graph->nvtxs; xadj = graph->xadj; adjncy = graph->adjncy; vwgt = graph->vwgt; adjwgt = graph->adjwgt; where = graph->where; pwgts = graph->pwgts; /* We assume that this is properly initialized */ maxpwgt = idxwspacemalloc(ctrl, nparts); ndoms = idxwspacemalloc(ctrl, nparts); otherpmat = idxwspacemalloc(ctrl, nparts); ind = idxwspacemalloc(ctrl, nvtxs); pmat = ctrl->wspace.pmat; cand = (KeyValueType *)GKmalloc(nparts*sizeof(KeyValueType), "EliminateSubDomainEdges: cand"); cand2 = (KeyValueType *)GKmalloc(nparts*sizeof(KeyValueType), "EliminateSubDomainEdges: cand"); /* Compute the pmat matrix and ndoms */ ComputeSubDomainGraph(graph, nparts, pmat, ndoms); /* Compute the maximum allowed weight for each domain */ tvwgt = idxsum(nparts, pwgts); for (i=0; i<nparts; i++) maxpwgt[i] = 1.25*tpwgts[i]*tvwgt; /* Get into the loop eliminating subdomain connections */ for (;;) { total = idxsum(nparts, ndoms); avg = total/nparts; max = ndoms[idxamax(nparts, ndoms)]; /* printf("Adjacent Subdomain Stats: Total: %3d, Max: %3d, Avg: %3d [%5d]\n", total, max, avg, idxsum(nparts*nparts, pmat)); */ if (max < 1.4*avg) break; me = idxamax(nparts, ndoms); mypmat = pmat + me*nparts; totalout = idxsum(nparts, mypmat); /*printf("Me: %d, TotalOut: %d,\n", me, totalout);*/ /* Sort the connections according to their cut */ for (ncand2=0, i=0; i<nparts; i++) { if (mypmat[i] > 0) { cand2[ncand2].key = mypmat[i]; cand2[ncand2++].val = i; } } ikeysort(ncand2, cand2); move = 0; for (min=0; min<ncand2; min++) { if (cand2[min].key > totalout/(2*ndoms[me])) break; other = cand2[min].val; /*printf("\tMinOut: %d to %d\n", mypmat[other], other);*/ idxset(nparts, 0, otherpmat); /* Go and find the vertices in 'other' that are connected in 'me' */ for (nind=0, i=0; i<nvtxs; i++) { if (where[i] == other) { for (j=xadj[i]; j<xadj[i+1]; j++) { if (where[adjncy[j]] == me) { ind[nind++] = i; break; } } } } /* Go and construct the otherpmat to see where these nind vertices are connected to */ for (cpwgt=0, ii=0; ii<nind; ii++) { i = ind[ii]; cpwgt += vwgt[i]; for (j=xadj[i]; j<xadj[i+1]; j++) otherpmat[where[adjncy[j]]] += adjwgt[j]; } otherpmat[other] = 0; for (ncand=0, i=0; i<nparts; i++) { if (otherpmat[i] > 0) { cand[ncand].key = -otherpmat[i]; cand[ncand++].val = i; } } ikeysort(ncand, cand); /* * Go through and the select the first domain that is common with 'me', and * does not increase the ndoms[target] higher than my ndoms, subject to the * maxpwgt constraint. Traversal is done from the mostly connected to the least. */ target = target2 = -1; for (i=0; i<ncand; i++) { k = cand[i].val; if (mypmat[k] > 0) { if (pwgts[k] + cpwgt > maxpwgt[k]) /* Check if balance will go off */ continue; for (j=0; j<nparts; j++) { if (otherpmat[j] > 0 && ndoms[j] >= ndoms[me]-1 && pmat[nparts*j+k] == 0) break; } if (j == nparts) { /* No bad second level effects */ for (nadd=0, j=0; j<nparts; j++) { if (otherpmat[j] > 0 && pmat[nparts*k+j] == 0) nadd++; } /*printf("\t\tto=%d, nadd=%d, %d\n", k, nadd, ndoms[k]);*/ if (target2 == -1 && ndoms[k]+nadd < ndoms[me]) { target2 = k; } if (nadd == 0) { target = k; break; } } } } if (target == -1 && target2 != -1) target = target2; if (target == -1) { /* printf("\t\tCould not make the move\n");*/ continue; } /*printf("\t\tMoving to %d\n", target);*/ /* Update the partition weights */ INC_DEC(pwgts[target], pwgts[other], cpwgt); MoveGroupMConn(ctrl, graph, ndoms, pmat, nparts, target, nind, ind); move = 1; break; } if (move == 0) break; } idxwspacefree(ctrl, nparts); idxwspacefree(ctrl, nparts); idxwspacefree(ctrl, nparts); idxwspacefree(ctrl, nvtxs); GKfree(&cand, &cand2, LTERM); }
/************************************************************************* * This function compresses a graph by merging identical vertices * The compression should lead to at least 10% reduction. **************************************************************************/ void CompressGraph(CtrlType *ctrl, GraphType *graph, int nvtxs, idxtype *xadj, idxtype *adjncy, idxtype *cptr, idxtype *cind) { int i, ii, iii, j, jj, k, l, cnvtxs, cnedges; idxtype *cxadj, *cadjncy, *cvwgt, *mark, *map; KeyValueType *keys; mark = idxsmalloc(nvtxs, -1, "CompressGraph: mark"); map = idxsmalloc(nvtxs, -1, "CompressGraph: map"); keys = (KeyValueType *)GKmalloc(nvtxs*sizeof(KeyValueType), "CompressGraph: keys"); /* Compute a key for each adjacency list */ for (i=0; i<nvtxs; i++) { k = 0; for (j=xadj[i]; j<xadj[i+1]; j++) k += adjncy[j]; keys[i].key = k+i; /* Add the diagonal entry as well */ keys[i].val = i; } ikeysort(nvtxs, keys); l = cptr[0] = 0; for (cnvtxs=i=0; i<nvtxs; i++) { ii = keys[i].val; if (map[ii] == -1) { mark[ii] = i; /* Add the diagonal entry */ for (j=xadj[ii]; j<xadj[ii+1]; j++) mark[adjncy[j]] = i; cind[l++] = ii; map[ii] = cnvtxs; for (j=i+1; j<nvtxs; j++) { iii = keys[j].val; if (keys[i].key != keys[j].key || xadj[ii+1]-xadj[ii] != xadj[iii+1]-xadj[iii]) break; /* Break if keys or degrees are different */ if (map[iii] == -1) { /* Do a comparison if iii has not been mapped */ for (jj=xadj[iii]; jj<xadj[iii+1]; jj++) { if (mark[adjncy[jj]] != i) break; } if (jj == xadj[iii+1]) { /* Identical adjacency structure */ map[iii] = cnvtxs; cind[l++] = iii; } } } cptr[++cnvtxs] = l; } } /* printf("Original: %6d, Compressed: %6d\n", nvtxs, cnvtxs); */ InitGraph(graph); if (cnvtxs >= COMPRESSION_FRACTION*nvtxs) { graph->nvtxs = nvtxs; graph->nedges = xadj[nvtxs]; graph->ncon = 1; graph->xadj = xadj; graph->adjncy = adjncy; graph->gdata = idxmalloc(3*nvtxs+graph->nedges, "CompressGraph: gdata"); graph->vwgt = graph->gdata; graph->adjwgtsum = graph->gdata+nvtxs; graph->cmap = graph->gdata+2*nvtxs; graph->adjwgt = graph->gdata+3*nvtxs; idxset(nvtxs, 1, graph->vwgt); idxset(graph->nedges, 1, graph->adjwgt); for (i=0; i<nvtxs; i++) graph->adjwgtsum[i] = xadj[i+1]-xadj[i]; graph->label = idxmalloc(nvtxs, "CompressGraph: label"); for (i=0; i<nvtxs; i++) graph->label[i] = i; } else { /* Ok, form the compressed graph */ cnedges = 0; for (i=0; i<cnvtxs; i++) { ii = cind[cptr[i]]; cnedges += xadj[ii+1]-xadj[ii]; } /* Allocate memory for the compressed graph*/ graph->gdata = idxmalloc(4*cnvtxs+1 + 2*cnedges, "CompressGraph: gdata"); cxadj = graph->xadj = graph->gdata; cvwgt = graph->vwgt = graph->gdata + cnvtxs+1; graph->adjwgtsum = graph->gdata + 2*cnvtxs+1; graph->cmap = graph->gdata + 3*cnvtxs+1; cadjncy = graph->adjncy = graph->gdata + 4*cnvtxs+1; graph->adjwgt = graph->gdata + 4*cnvtxs+1 + cnedges; /* Now go and compress the graph */ idxset(nvtxs, -1, mark); l = cxadj[0] = 0; for (i=0; i<cnvtxs; i++) { cvwgt[i] = cptr[i+1]-cptr[i]; mark[i] = i; /* Remove any dioganal entries in the compressed graph */ for (j=cptr[i]; j<cptr[i+1]; j++) { ii = cind[j]; for (jj=xadj[ii]; jj<xadj[ii+1]; jj++) { k = map[adjncy[jj]]; if (mark[k] != i) cadjncy[l++] = k; mark[k] = i; } } cxadj[i+1] = l; } graph->nvtxs = cnvtxs; graph->nedges = l; graph->ncon = 1; idxset(graph->nedges, 1, graph->adjwgt); for (i=0; i<cnvtxs; i++) graph->adjwgtsum[i] = cxadj[i+1]-cxadj[i]; graph->label = idxmalloc(cnvtxs, "CompressGraph: label"); for (i=0; i<cnvtxs; i++) graph->label[i] = i; } GKfree(&keys, &map, &mark, LTERM); }
/************************************************************************* * This function performs an edge-based FM refinement **************************************************************************/ void Mc_Serial_Balance2Way(GraphType *graph, float *tpwgts, float lbfactor) { int i, ii, j, k, kwgt, nvtxs, ncon, nbnd, nswaps, from, to, limit, tmp, cnum; idxtype *xadj, *adjncy, *adjwgt, *where, *id, *ed, *bndptr, *bndind; idxtype *moved, *swaps, *qnum; float *nvwgt, *npwgts, mindiff[MAXNCON], origbal, minbal, newbal; FPQueueType parts[MAXNCON][2]; int higain, oldgain, mincut, newcut, mincutorder; int qsizes[MAXNCON][2]; KeyValueType *cand; nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; nvwgt = graph->nvwgt; adjncy = graph->adjncy; adjwgt = graph->adjwgt; where = graph->where; id = graph->sendind; ed = graph->recvind; npwgts = graph->gnpwgts; bndptr = graph->sendptr; bndind = graph->recvptr; moved = idxmalloc(nvtxs, "moved"); swaps = idxmalloc(nvtxs, "swaps"); qnum = idxmalloc(nvtxs, "qnum"); cand = (KeyValueType *)GKmalloc(nvtxs*sizeof(KeyValueType), "cand"); limit = amin(amax(0.01*nvtxs, 15), 100); /* Initialize the queues */ for (i=0; i<ncon; i++) { FPQueueInit(&parts[i][0], nvtxs); FPQueueInit(&parts[i][1], nvtxs); qsizes[i][0] = qsizes[i][1] = 0; } for (i=0; i<nvtxs; i++) { qnum[i] = samax(ncon, nvwgt+i*ncon); qsizes[qnum[i]][where[i]]++; } for (from=0; from<2; from++) { for (j=0; j<ncon; j++) { if (qsizes[j][from] == 0) { for (i=0; i<nvtxs; i++) { if (where[i] != from) continue; k = samax2(ncon, nvwgt+i*ncon); if (k == j && qsizes[qnum[i]][from] > qsizes[j][from] && nvwgt[i*ncon+qnum[i]] < 1.3*nvwgt[i*ncon+j]) { qsizes[qnum[i]][from]--; qsizes[j][from]++; qnum[i] = j; } } } } } for (i=0; i<ncon; i++) mindiff[i] = fabs(tpwgts[i]-npwgts[i]); minbal = origbal = Serial_Compute2WayHLoadImbalance(ncon, npwgts, tpwgts); newcut = mincut = graph->mincut; mincutorder = -1; idxset(nvtxs, -1, moved); /* Insert all nodes in the priority queues */ nbnd = graph->gnvtxs; for (i=0; i<nvtxs; i++) { cand[i].key = id[i]-ed[i]; cand[i].val = i; } ikeysort(nvtxs, cand); for (ii=0; ii<nvtxs; ii++) { i = cand[ii].val; FPQueueInsert(&parts[qnum[i]][where[i]], i, (float)(ed[i]-id[i])); } for (nswaps=0; nswaps<nvtxs; nswaps++) { if (minbal < lbfactor) break; Serial_SelectQueue(ncon, npwgts, tpwgts, &from, &cnum, parts); to = (from+1)%2; if (from == -1 || (higain = FPQueueGetMax(&parts[cnum][from])) == -1) break; saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+from*ncon, 1); newcut -= (ed[higain]-id[higain]); newbal = Serial_Compute2WayHLoadImbalance(ncon, npwgts, tpwgts); if (newbal < minbal || (newbal == minbal && (newcut < mincut || (newcut == mincut && Serial_BetterBalance(ncon, npwgts, tpwgts, mindiff))))) { mincut = newcut; minbal = newbal; mincutorder = nswaps; for (i=0; i<ncon; i++) mindiff[i] = fabs(tpwgts[i]-npwgts[i]); } else if (nswaps-mincutorder > limit) { /* We hit the limit, undo last move */ newcut += (ed[higain]-id[higain]); saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+from*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); break; } where[higain] = to; moved[higain] = nswaps; swaps[nswaps] = higain; /************************************************************** * Update the id[i]/ed[i] values of the affected nodes ***************************************************************/ SWAP(id[higain], ed[higain], tmp); if (ed[higain] == 0 && bndptr[higain] != -1 && xadj[higain] < xadj[higain+1]) BNDDelete(nbnd, bndind, bndptr, higain); if (ed[higain] > 0 && bndptr[higain] == -1) BNDInsert(nbnd, bndind, bndptr, higain); for (j=xadj[higain]; j<xadj[higain+1]; j++) { k = adjncy[j]; oldgain = ed[k]-id[k]; kwgt = (to == where[k] ? adjwgt[j] : -adjwgt[j]); INC_DEC(id[k], ed[k], kwgt); /* Update the queue position */ if (moved[k] == -1) FPQueueUpdate(&parts[qnum[k]][where[k]], k, (float)(oldgain), (float)(ed[k]-id[k])); /* Update its boundary information */ if (ed[k] == 0 && bndptr[k] != -1) BNDDelete(nbnd, bndind, bndptr, k); else if (ed[k] > 0 && bndptr[k] == -1) BNDInsert(nbnd, bndind, bndptr, k); } } /**************************************************************** * Roll back computations *****************************************************************/ for (nswaps--; nswaps>mincutorder; nswaps--) { higain = swaps[nswaps]; to = where[higain] = (where[higain]+1)%2; SWAP(id[higain], ed[higain], tmp); if (ed[higain] == 0 && bndptr[higain] != -1 && xadj[higain] < xadj[higain+1]) BNDDelete(nbnd, bndind, bndptr, higain); else if (ed[higain] > 0 && bndptr[higain] == -1) BNDInsert(nbnd, bndind, bndptr, higain); saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+((to+1)%2)*ncon, 1); for (j=xadj[higain]; j<xadj[higain+1]; j++) { k = adjncy[j]; kwgt = (to == where[k] ? adjwgt[j] : -adjwgt[j]); INC_DEC(id[k], ed[k], kwgt); if (bndptr[k] != -1 && ed[k] == 0) BNDDelete(nbnd, bndind, bndptr, k); if (bndptr[k] == -1 && ed[k] > 0) BNDInsert(nbnd, bndind, bndptr, k); } } graph->mincut = mincut; graph->gnvtxs = nbnd; for (i=0; i<ncon; i++) { FPQueueFree(&parts[i][0]); FPQueueFree(&parts[i][1]); } GKfree((void **)&cand, (void **)&qnum, (void **)&moved, (void **)&swaps, LTERM); return; }
/************************************************************************* * This function performs an edge-based FM refinement **************************************************************************/ void Mc_Serial_FM_2WayRefine(GraphType *graph, float *tpwgts, int npasses) { int i, ii, j, k; int kwgt, nvtxs, ncon, nbnd, nswaps, from, to, pass, limit, tmp, cnum; idxtype *xadj, *adjncy, *adjwgt, *where, *id, *ed, *bndptr, *bndind; idxtype *moved, *swaps, *qnum; float *nvwgt, *npwgts, mindiff[MAXNCON], origbal, minbal, newbal; FPQueueType parts[MAXNCON][2]; int higain, oldgain, mincut, initcut, newcut, mincutorder; float rtpwgts[MAXNCON*2]; KeyValueType *cand; int mype; MPI_Comm_rank(MPI_COMM_WORLD, &mype); nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; nvwgt = graph->nvwgt; adjncy = graph->adjncy; adjwgt = graph->adjwgt; where = graph->where; id = graph->sendind; ed = graph->recvind; npwgts = graph->gnpwgts; bndptr = graph->sendptr; bndind = graph->recvptr; moved = idxmalloc(nvtxs, "moved"); swaps = idxmalloc(nvtxs, "swaps"); qnum = idxmalloc(nvtxs, "qnum"); cand = (KeyValueType *)GKmalloc(nvtxs*sizeof(KeyValueType), "cand"); limit = amin(amax(0.01*nvtxs, 25), 150); /* Initialize the queues */ for (i=0; i<ncon; i++) { FPQueueInit(&parts[i][0], nvtxs); FPQueueInit(&parts[i][1], nvtxs); } for (i=0; i<nvtxs; i++) qnum[i] = samax(ncon, nvwgt+i*ncon); origbal = Serial_Compute2WayHLoadImbalance(ncon, npwgts, tpwgts); for (i=0; i<ncon; i++) { rtpwgts[i] = origbal*tpwgts[i]; rtpwgts[ncon+i] = origbal*tpwgts[ncon+i]; } idxset(nvtxs, -1, moved); for (pass=0; pass<npasses; pass++) { /* Do a number of passes */ for (i=0; i<ncon; i++) { FPQueueReset(&parts[i][0]); FPQueueReset(&parts[i][1]); } mincutorder = -1; newcut = mincut = initcut = graph->mincut; for (i=0; i<ncon; i++) mindiff[i] = fabs(tpwgts[i]-npwgts[i]); minbal = Serial_Compute2WayHLoadImbalance(ncon, npwgts, tpwgts); /* Insert boundary nodes in the priority queues */ nbnd = graph->gnvtxs; for (i=0; i<nbnd; i++) { cand[i].key = id[i]-ed[i]; cand[i].val = i; } ikeysort(nbnd, cand); for (ii=0; ii<nbnd; ii++) { i = bndind[cand[ii].val]; FPQueueInsert(&parts[qnum[i]][where[i]], i, (float)(ed[i]-id[i])); } for (nswaps=0; nswaps<nvtxs; nswaps++) { Serial_SelectQueue(ncon, npwgts, rtpwgts, &from, &cnum, parts); to = (from+1)%2; if (from == -1 || (higain = FPQueueGetMax(&parts[cnum][from])) == -1) break; saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+from*ncon, 1); newcut -= (ed[higain]-id[higain]); newbal = Serial_Compute2WayHLoadImbalance(ncon, npwgts, tpwgts); if ((newcut < mincut && newbal-origbal <= .00001) || (newcut == mincut && (newbal < minbal || (newbal == minbal && Serial_BetterBalance(ncon, npwgts, tpwgts, mindiff))))) { mincut = newcut; minbal = newbal; mincutorder = nswaps; for (i=0; i<ncon; i++) mindiff[i] = fabs(tpwgts[i]-npwgts[i]); } else if (nswaps-mincutorder > limit) { /* We hit the limit, undo last move */ newcut += (ed[higain]-id[higain]); saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+from*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); break; } where[higain] = to; moved[higain] = nswaps; swaps[nswaps] = higain; /************************************************************** * Update the id[i]/ed[i] values of the affected nodes ***************************************************************/ SWAP(id[higain], ed[higain], tmp); if (ed[higain] == 0 && xadj[higain] < xadj[higain+1]) BNDDelete(nbnd, bndind, bndptr, higain); for (j=xadj[higain]; j<xadj[higain+1]; j++) { k = adjncy[j]; oldgain = ed[k]-id[k]; kwgt = (to == where[k] ? adjwgt[j] : -adjwgt[j]); INC_DEC(id[k], ed[k], kwgt); /* Update its boundary information and queue position */ if (bndptr[k] != -1) { /* If k was a boundary vertex */ if (ed[k] == 0) { /* Not a boundary vertex any more */ BNDDelete(nbnd, bndind, bndptr, k); if (moved[k] == -1) /* Remove it if in the queues */ FPQueueDelete(&parts[qnum[k]][where[k]], k); } else { /* If it has not been moved, update its position in the queue */ if (moved[k] == -1) FPQueueUpdate(&parts[qnum[k]][where[k]], k, (float)oldgain, (float)(ed[k]-id[k])); } } else { if (ed[k] > 0) { /* It will now become a boundary vertex */ BNDInsert(nbnd, bndind, bndptr, k); if (moved[k] == -1) FPQueueInsert(&parts[qnum[k]][where[k]], k, (float)(ed[k]-id[k])); } } } } /**************************************************************** * Roll back computations *****************************************************************/ for (i=0; i<nswaps; i++) moved[swaps[i]] = -1; /* reset moved array */ for (nswaps--; nswaps>mincutorder; nswaps--) { higain = swaps[nswaps]; to = where[higain] = (where[higain]+1)%2; SWAP(id[higain], ed[higain], tmp); if (ed[higain] == 0 && bndptr[higain] != -1 && xadj[higain] < xadj[higain+1]) BNDDelete(nbnd, bndind, bndptr, higain); else if (ed[higain] > 0 && bndptr[higain] == -1) BNDInsert(nbnd, bndind, bndptr, higain); saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+((to+1)%2)*ncon, 1); for (j=xadj[higain]; j<xadj[higain+1]; j++) { k = adjncy[j]; kwgt = (to == where[k] ? adjwgt[j] : -adjwgt[j]); INC_DEC(id[k], ed[k], kwgt); if (bndptr[k] != -1 && ed[k] == 0) BNDDelete(nbnd, bndind, bndptr, k); if (bndptr[k] == -1 && ed[k] > 0) BNDInsert(nbnd, bndind, bndptr, k); } } graph->mincut = mincut; graph->gnvtxs = nbnd; if (mincutorder == -1 || mincut == initcut) break; } for (i=0; i<ncon; i++) { FPQueueFree(&parts[i][0]); FPQueueFree(&parts[i][1]); } GKfree((void **)&cand, (void **)&qnum, (void **)&moved, (void **)&swaps, LTERM); return; }
/************************************************************** * This subroutine remaps a partitioning on a single processor **************************************************************/ void SerialRemap(GraphType *graph, int nparts, idxtype *base, idxtype *scratch, idxtype *remap, float *tpwgts) { int i, ii, j, k; int nvtxs, nmapped, max_mult; int from, to, current_from, smallcount, bigcount; KeyValueType *flowto, *bestflow; KeyKeyValueType *sortvtx; idxtype *vsize, *htable, *map, *rowmap; nvtxs = graph->nvtxs; vsize = graph->vsize; max_mult = amin(MAX_NPARTS_MULTIPLIER, nparts); sortvtx = (KeyKeyValueType *)GKmalloc(nvtxs*sizeof(KeyKeyValueType), "sortvtx"); flowto = (KeyValueType *)GKmalloc((nparts*max_mult+nparts)*sizeof(KeyValueType), "flowto"); bestflow = flowto+nparts; map = htable = idxsmalloc(nparts*2, -1, "htable"); rowmap = map+nparts; for (i=0; i<nvtxs; i++) { sortvtx[i].key1 = base[i]; sortvtx[i].key2 = vsize[i]; sortvtx[i].val = i; } qsort((void *)sortvtx, (size_t)nvtxs, (size_t)sizeof(KeyKeyValueType), SSMIncKeyCmp); for (j=0; j<nparts; j++) { flowto[j].key = 0; flowto[j].val = j; } /* this step has nparts*nparts*log(nparts) computational complexity */ bigcount = smallcount = current_from = 0; for (ii=0; ii<nvtxs; ii++) { i = sortvtx[ii].val; from = base[i]; to = scratch[i]; if (from > current_from) { /* reset the hash table */ for (j=0; j<smallcount; j++) htable[flowto[j].val] = -1; ASSERTS(idxsum(nparts, htable) == -nparts); ikeysort(smallcount, flowto); for (j=0; j<amin(smallcount, max_mult); j++, bigcount++) { bestflow[bigcount].key = flowto[j].key; bestflow[bigcount].val = current_from*nparts+flowto[j].val; } smallcount = 0; current_from = from; } if (htable[to] == -1) { htable[to] = smallcount; flowto[smallcount].key = -vsize[i]; flowto[smallcount].val = to; smallcount++; } else { flowto[htable[to]].key += -vsize[i]; } } /* reset the hash table */ for (j=0; j<smallcount; j++) htable[flowto[j].val] = -1; ASSERTS(idxsum(nparts, htable) == -nparts); ikeysort(smallcount, flowto); for (j=0; j<amin(smallcount, max_mult); j++, bigcount++) { bestflow[bigcount].key = flowto[j].key; bestflow[bigcount].val = current_from*nparts+flowto[j].val; } ikeysort(bigcount, bestflow); ASSERTS(idxsum(nparts, map) == -nparts); ASSERTS(idxsum(nparts, rowmap) == -nparts); nmapped = 0; /* now make as many assignments as possible */ for (ii=0; ii<bigcount; ii++) { i = bestflow[ii].val; j = i % nparts; /* to */ k = i / nparts; /* from */ if (map[j] == -1 && rowmap[k] == -1 && SimilarTpwgts(tpwgts, graph->ncon, j, k)) { map[j] = k; rowmap[k] = j; nmapped++; } if (nmapped == nparts) break; } /* remap the rest */ /* it may help try remapping to the same label first */ if (nmapped < nparts) { for (j=0; j<nparts && nmapped<nparts; j++) { if (map[j] == -1) { for (ii=0; ii<nparts; ii++) { i = (j+ii) % nparts; if (rowmap[i] == -1 && SimilarTpwgts(tpwgts, graph->ncon, i, j)) { map[j] = i; rowmap[i] = j; nmapped++; break; } } } } } /* check to see if remapping fails (due to dis-similar tpwgts) */ /* if remapping fails, revert to original mapping */ if (nmapped < nparts) for (i=0; i<nparts; i++) map[i] = i; for (i=0; i<nvtxs; i++) remap[i] = map[remap[i]]; GKfree((void **)&sortvtx, (void **)&flowto, (void **)&htable, LTERM); }
/************************************************************************* * 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; }
/************************************************************************* * This function balances two partitions by moving the highest gain * (including negative gain) vertices to the other domain. * It is used only when tha unbalance is due to non contigous * subdomains. That is, the are no boundary vertices. * It moves vertices from the domain that is overweight to the one that * is underweight. **************************************************************************/ void Mc_Serial_Init2WayBalance(GraphType *graph, float *tpwgts) { int i, ii, j, k; int kwgt, nvtxs, nbnd, ncon, nswaps, from, to, cnum, tmp; idxtype *xadj, *adjncy, *adjwgt, *where, *id, *ed, *bndptr, *bndind; idxtype *qnum; float *nvwgt, *npwgts; FPQueueType parts[MAXNCON][2]; int higain, oldgain, mincut; KeyValueType *cand; nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; adjncy = graph->adjncy; nvwgt = graph->nvwgt; adjwgt = graph->adjwgt; where = graph->where; id = graph->sendind; ed = graph->recvind; npwgts = graph->gnpwgts; bndptr = graph->sendptr; bndind = graph->recvptr; qnum = idxmalloc(nvtxs, "qnum"); cand = (KeyValueType *)GKmalloc(nvtxs*sizeof(KeyValueType), "cand"); /* This is called for initial partitioning so we know from where to pick nodes */ from = 1; to = (from+1)%2; for (i=0; i<ncon; i++) { FPQueueInit(&parts[i][0], nvtxs); FPQueueInit(&parts[i][1], nvtxs); } /* Compute the queues in which each vertex will be assigned to */ for (i=0; i<nvtxs; i++) qnum[i] = samax(ncon, nvwgt+i*ncon); for (i=0; i<nvtxs; i++) { cand[i].key = id[i]-ed[i]; cand[i].val = i; } ikeysort(nvtxs, cand); /* Insert the nodes of the proper partition in the appropriate priority queue */ for (ii=0; ii<nvtxs; ii++) { i = cand[ii].val; if (where[i] == from) { if (ed[i] > 0) FPQueueInsert(&parts[qnum[i]][0], i, (float)(ed[i]-id[i])); else FPQueueInsert(&parts[qnum[i]][1], i, (float)(ed[i]-id[i])); } } mincut = graph->mincut; nbnd = graph->gnvtxs; for (nswaps=0; nswaps<nvtxs; nswaps++) { if (Serial_AreAnyVwgtsBelow(ncon, 1.0, npwgts+from*ncon, 0.0, nvwgt, tpwgts+from*ncon)) break; if ((cnum = Serial_SelectQueueOneWay(ncon, npwgts, tpwgts, from, parts)) == -1) break; if ((higain = FPQueueGetMax(&parts[cnum][0])) == -1) higain = FPQueueGetMax(&parts[cnum][1]); mincut -= (ed[higain]-id[higain]); saxpy2(ncon, 1.0, nvwgt+higain*ncon, 1, npwgts+to*ncon, 1); saxpy2(ncon, -1.0, nvwgt+higain*ncon, 1, npwgts+from*ncon, 1); where[higain] = to; /************************************************************** * Update the id[i]/ed[i] values of the affected nodes ***************************************************************/ SWAP(id[higain], ed[higain], tmp); if (ed[higain] == 0 && bndptr[higain] != -1 && xadj[higain] < xadj[higain+1]) BNDDelete(nbnd, bndind, bndptr, higain); if (ed[higain] > 0 && bndptr[higain] == -1) BNDInsert(nbnd, bndind, bndptr, higain); for (j=xadj[higain]; j<xadj[higain+1]; j++) { k = adjncy[j]; oldgain = ed[k]-id[k]; kwgt = (to == where[k] ? adjwgt[j] : -adjwgt[j]); INC_DEC(id[k], ed[k], kwgt); /* Update the queue position */ if (where[k] == from) { if (ed[k] > 0 && bndptr[k] == -1) { /* It moves in boundary */ FPQueueDelete(&parts[qnum[k]][1], k); FPQueueInsert(&parts[qnum[k]][0], k, (float)(ed[k]-id[k])); } else { /* It must be in the boundary already */ FPQueueUpdate(&parts[qnum[k]][0], k, (float)(oldgain), (float)(ed[k]-id[k])); } } /* Update its boundary information */ if (ed[k] == 0 && bndptr[k] != -1) BNDDelete(nbnd, bndind, bndptr, k); else if (ed[k] > 0 && bndptr[k] == -1) BNDInsert(nbnd, bndind, bndptr, k); } } graph->mincut = mincut; graph->gnvtxs = nbnd; for (i=0; i<ncon; i++) { FPQueueFree(&parts[i][0]); FPQueueFree(&parts[i][1]); } GKfree((void **)&cand, (void **)&qnum, LTERM); }
/************************************************************************* * This function computes the assignment using the the objective the * minimization of the total volume of data that needs to move **************************************************************************/ void ParallelTotalVReMap(CtrlType *ctrl, idxtype *lpwgts, idxtype *map, WorkSpaceType *wspace, int npasses, int ncon) { int i, ii, j, k, nparts, mype; int pass, maxipwgt, nmapped, oldwgt, newwgt, done; idxtype *rowmap, *mylpwgts; KeyValueType *recv, send; int nsaved, gnsaved; mype = ctrl->mype; nparts = ctrl->nparts; recv = (KeyValueType *)GKmalloc(sizeof(KeyValueType)*nparts, "remap: recv"); mylpwgts = idxmalloc(nparts, "mylpwgts"); done = nmapped = 0; idxset(nparts, -1, map); rowmap = idxset(nparts, -1, wspace->pv3); idxcopy(nparts, lpwgts, mylpwgts); for (pass=0; pass<npasses; pass++) { maxipwgt = idxamax(nparts, mylpwgts); if (mylpwgts[maxipwgt] > 0 && !done) { send.key = -mylpwgts[maxipwgt]; send.val = mype*nparts+maxipwgt; } else { send.key = 0; send.val = -1; } /* each processor sends its selection */ MPI_Allgather((void *)&send, 2, IDX_DATATYPE, (void *)recv, 2, IDX_DATATYPE, ctrl->comm); ikeysort(nparts, recv); if (recv[0].key == 0) break; /* now make as many assignments as possible */ for (ii=0; ii<nparts; ii++) { i = recv[ii].val; if (i == -1) continue; j = i % nparts; k = i / nparts; if (map[j] == -1 && rowmap[k] == -1 && SimilarTpwgts(ctrl->tpwgts, ncon, j, k)) { map[j] = k; rowmap[k] = j; nmapped++; mylpwgts[j] = 0; if (mype == k) done = 1; } if (nmapped == nparts) break; } if (nmapped == nparts) break; } /* Map unmapped partitions */ if (nmapped < nparts) { for (i=j=0; j<nparts && nmapped<nparts; j++) { if (map[j] == -1) { for (; i<nparts; i++) { if (rowmap[i] == -1 && SimilarTpwgts(ctrl->tpwgts, ncon, i, j)) { map[j] = i; rowmap[i] = j; nmapped++; break; } } } } } /* check to see if remapping fails (due to dis-similar tpwgts) */ /* if remapping fails, revert to original mapping */ if (nmapped < nparts) { for (i=0; i<nparts; i++) map[i] = i; IFSET(ctrl->dbglvl, DBG_REMAP, rprintf(ctrl, "Savings from parallel remapping: %0\n")); } else { /* check for a savings */ oldwgt = lpwgts[mype]; newwgt = lpwgts[rowmap[mype]]; nsaved = newwgt - oldwgt; gnsaved = GlobalSESum(ctrl, nsaved); /* undo everything if we don't see a savings */ if (gnsaved <= 0) { for (i=0; i<nparts; i++) map[i] = i; } IFSET(ctrl->dbglvl, DBG_REMAP, rprintf(ctrl, "Savings from parallel remapping: %d\n", amax(0,gnsaved))); } GKfree((void **)&recv, (void **)&mylpwgts, LTERM); }
/************************************************************************* * This function tests the repeated shmem_put **************************************************************************/ void SetUp(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace) { int i, j, k, islocal, penum, gnvtxs, nvtxs, nlocal, firstvtx, lastvtx, nsend, nrecv, nnbrs, nadj; int npes=ctrl->npes, mype=ctrl->mype; idxtype *vtxdist, *xadj, *adjncy; idxtype *peind, *recvptr, *recvind, *sendptr, *sendind; idxtype *receive, *pemap, *imap, *lperm; idxtype *pexadj, *peadjncy, *peadjloc, *startsind; KeyValueType *recvrequests, *sendrequests, *adjpairs; IFSET(ctrl->dbglvl, DBG_TIME, MPI_Barrier(ctrl->comm)); IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->SetupTmr)); gnvtxs = graph->gnvtxs; nvtxs = graph->nvtxs; vtxdist = graph->vtxdist; xadj = graph->xadj; adjncy = graph->adjncy; firstvtx = vtxdist[mype]; lastvtx = vtxdist[mype+1]; pemap = wspace->pv1; idxset(npes, -1, pemap); lperm = graph->lperm = idxmalloc(nvtxs, "SetUp: graph->lperm"); for (i=0; i<nvtxs; i++) lperm[i] = i; /************************************************************* * Determine what you need to receive *************************************************************/ receive = wspace->indices; /* Use the large global received array for now */ adjpairs = wspace->pairs; for (nlocal = nadj = i = 0; i<nvtxs; i++) { islocal = 1; for (j=xadj[i]; j<xadj[i+1]; j++) { k = adjncy[j]; if (k >= firstvtx && k < lastvtx) { adjncy[j] = k-firstvtx; continue; /* local vertex */ } adjpairs[nadj].key = k; adjpairs[nadj++].val = j; islocal = 0; } if (islocal) { lperm[i] = lperm[nlocal]; lperm[nlocal++] = i; } } /* Take care the received part now */ ikeysort(nadj, adjpairs); adjpairs[nadj].key = gnvtxs+1; /* Boundary condition */ for (nrecv=i=0; i<nadj; i++) { adjncy[adjpairs[i].val] = nvtxs+nrecv; if (adjpairs[i].key != adjpairs[i+1].key) receive[nrecv++] = adjpairs[i].key; } /* Allocate space for the setup info attached to this level of the graph */ peind = graph->peind = idxmalloc(npes, "SetUp: peind"); recvptr = graph->recvptr = idxmalloc(npes+1, "SetUp: recvptr"); recvind = graph->recvind = idxmalloc(nrecv, "SetUp: recvind"); /* Take care of the received portion */ idxcopy(nrecv, receive, recvind); /* Copy the vertices to be received into recvind */ i = nnbrs = recvptr[0] = 0; for (penum=0; penum<npes; penum++) { for (j=i; j<nrecv; j++) { if (recvind[j] >= vtxdist[penum+1]) break; } if (j > i) { peind[nnbrs] = penum; recvptr[++nnbrs] = j; i = j; } } /************************************************************* * Determine what you need to send *************************************************************/ /* Tell the other processors what they need to send you */ recvrequests = wspace->pepairs1; sendrequests = wspace->pepairs2; for (i=0; i<npes; i++) recvrequests[i].key = 0; for (i=0; i<nnbrs; i++) { recvrequests[peind[i]].key = recvptr[i+1]-recvptr[i]; recvrequests[peind[i]].val = nvtxs+recvptr[i]; } MPI_Alltoall((void *)recvrequests, 2, IDX_DATATYPE, (void *)sendrequests, 2, IDX_DATATYPE, ctrl->comm); sendptr = graph->sendptr = idxmalloc(npes+1, "SetUp: sendptr"); startsind = wspace->pv2; for (j=i=0; i<npes; i++) { if (sendrequests[i].key > 0) { sendptr[j] = sendrequests[i].key; startsind[j] = sendrequests[i].val; j++; } } ASSERT(ctrl, nnbrs == j); MAKECSR(i, j, sendptr); nsend = sendptr[nnbrs]; sendind = graph->sendind = idxmalloc(nsend, "SetUp: sendind"); /* Issue the receives for sendind */ for (i=0; i<nnbrs; i++) { MPI_Irecv((void *)(sendind+sendptr[i]), sendptr[i+1]-sendptr[i], IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->rreq+i); } /* Issue the sends. My recvind[penum] becomes penum's sendind[mype] */ for (i=0; i<nnbrs; i++) { MPI_Isend((void *)(recvind+recvptr[i]), recvptr[i+1]-recvptr[i], IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->sreq+i); } MPI_Waitall(nnbrs, ctrl->rreq, ctrl->statuses); MPI_Waitall(nnbrs, ctrl->sreq, ctrl->statuses); /* Create the peadjncy data structure for sparse boundary exchanges */ pexadj = graph->pexadj = idxsmalloc(nvtxs+1, 0, "SetUp: pexadj"); peadjncy = graph->peadjncy = idxmalloc(nsend, "SetUp: peadjncy"); peadjloc = graph->peadjloc = idxmalloc(nsend, "SetUp: peadjloc"); for (i=0; i<nsend; i++) { ASSERTP(ctrl, sendind[i] >= firstvtx && sendind[i] < lastvtx, (ctrl, "%d %d %d\n", sendind[i], firstvtx, lastvtx)); pexadj[sendind[i]-firstvtx]++; } MAKECSR(i, nvtxs, pexadj); for (i=0; i<nnbrs; i++) { for (j=sendptr[i]; j<sendptr[i+1]; j++) { k = pexadj[sendind[j]-firstvtx]++; peadjncy[k] = i; /* peind[i] is the actual PE number */ peadjloc[k] = startsind[i]++; } } ASSERT(ctrl, pexadj[nvtxs] == nsend); for (i=nvtxs; i>0; i--) pexadj[i] = pexadj[i-1]; pexadj[0] = 0; graph->nnbrs = nnbrs; graph->nrecv = nrecv; graph->nsend = nsend; graph->nlocal = nlocal; /* Create the inverse map from ladjncy to adjncy */ imap = graph->imap = idxmalloc(nvtxs+nrecv, "SetUp: imap"); for (i=0; i<nvtxs; i++) imap[i] = firstvtx+i; for (i=0; i<nrecv; i++) imap[nvtxs+i] = recvind[i]; /* Check if wspace->nlarge is large enough for nrecv and nsend */ if (wspace->nlarge < nrecv+nsend) { free(wspace->indices); free(wspace->pairs); wspace->nlarge = nrecv+nsend; wspace->indices = idxmalloc(wspace->nlarge, "SetUp: wspace->indices"); wspace->pairs = (KeyValueType *)GKmalloc(sizeof(KeyValueType)*wspace->nlarge, "SetUp: wspace->pairs"); } IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->SetupTmr)); #ifdef DEBUG_SETUPINFO rprintf(ctrl, "[%5d %5d] \tl:[%5d %5d] \ts:[%5d, %5d] \tr:[%5d, %5d]\n", GlobalSEMin(ctrl, nvtxs), GlobalSEMax(ctrl, nvtxs), GlobalSEMin(ctrl, nlocal), GlobalSEMax(ctrl, nlocal), GlobalSEMin(ctrl, nsend), GlobalSEMax(ctrl, nsend), GlobalSEMin(ctrl, nrecv), GlobalSEMax(ctrl, nrecv)); PrintSetUpInfo(ctrl, graph); #endif }
/************************************************************************* * This function converts a mesh into a dual graph **************************************************************************/ void ParMETIS_V3_Mesh2Dual(idxtype *elmdist, idxtype *eptr, idxtype *eind, int *numflag, int *ncommonnodes, idxtype **xadj, idxtype **adjncy, MPI_Comm *comm) { int i, j, jj, k, kk, m; int npes, mype, pe, count, mask, pass; int nelms, lnns, my_nns, node; int firstelm, firstnode, lnode, nrecv, nsend; int *scounts, *rcounts, *sdispl, *rdispl; idxtype *nodedist, *nmap, *auxarray; idxtype *gnptr, *gnind, *nptr, *nind, *myxadj, *myadjncy = NULL; idxtype *sbuffer, *rbuffer, *htable; KeyValueType *nodelist, *recvbuffer; idxtype ind[200], wgt[200]; int gmaxnode, gminnode; CtrlType ctrl; SetUpCtrl(&ctrl, -1, 0, *comm); npes = ctrl.npes; mype = ctrl.mype; nelms = elmdist[mype+1]-elmdist[mype]; if (*numflag == 1) ChangeNumberingMesh2(elmdist, eptr, eind, NULL, NULL, NULL, npes, mype, 1); mask = (1<<11)-1; /*****************************/ /* Determine number of nodes */ /*****************************/ gminnode = GlobalSEMin(&ctrl, eind[idxamin(eptr[nelms], eind)]); for (i=0; i<eptr[nelms]; i++) eind[i] -= gminnode; gmaxnode = GlobalSEMax(&ctrl, eind[idxamax(eptr[nelms], eind)]); /**************************/ /* Check for input errors */ /**************************/ ASSERTS(nelms > 0); /* construct node distribution array */ nodedist = idxsmalloc(npes+1, 0, "nodedist"); for (nodedist[0]=0, i=0,j=gmaxnode+1; i<npes; i++) { k = j/(npes-i); nodedist[i+1] = nodedist[i]+k; j -= k; } my_nns = nodedist[mype+1]-nodedist[mype]; firstnode = nodedist[mype]; nodelist = (KeyValueType *)GKmalloc(eptr[nelms]*sizeof(KeyValueType), "nodelist"); auxarray = idxmalloc(eptr[nelms], "auxarray"); htable = idxsmalloc(amax(my_nns, mask+1), -1, "htable"); scounts = imalloc(4*npes+2, "scounts"); rcounts = scounts+npes; sdispl = scounts+2*npes; rdispl = scounts+3*npes+1; /*********************************************/ /* first find a local numbering of the nodes */ /*********************************************/ for (i=0; i<nelms; i++) { for (j=eptr[i]; j<eptr[i+1]; j++) { nodelist[j].key = eind[j]; nodelist[j].val = j; auxarray[j] = i; /* remember the local element ID that uses this node */ } } ikeysort(eptr[nelms], nodelist); for (count=1, i=1; i<eptr[nelms]; i++) { if (nodelist[i].key > nodelist[i-1].key) count++; } lnns = count; nmap = idxmalloc(lnns, "nmap"); /* renumber the nodes of the elements array */ count = 1; nmap[0] = nodelist[0].key; eind[nodelist[0].val] = 0; nodelist[0].val = auxarray[nodelist[0].val]; /* Store the local element ID */ for (i=1; i<eptr[nelms]; i++) { if (nodelist[i].key > nodelist[i-1].key) { nmap[count] = nodelist[i].key; count++; } eind[nodelist[i].val] = count-1; nodelist[i].val = auxarray[nodelist[i].val]; /* Store the local element ID */ } MPI_Barrier(*comm); /**********************************************************/ /* perform comms necessary to construct node-element list */ /**********************************************************/ iset(npes, 0, scounts); for (pe=i=0; i<eptr[nelms]; i++) { while (nodelist[i].key >= nodedist[pe+1]) pe++; scounts[pe] += 2; } ASSERTS(pe < npes); MPI_Alltoall((void *)scounts, 1, MPI_INT, (void *)rcounts, 1, MPI_INT, *comm); icopy(npes, scounts, sdispl); MAKECSR(i, npes, sdispl); icopy(npes, rcounts, rdispl); MAKECSR(i, npes, rdispl); ASSERTS(sdispl[npes] == eptr[nelms]*2); nrecv = rdispl[npes]/2; recvbuffer = (KeyValueType *)GKmalloc(amax(1, nrecv)*sizeof(KeyValueType), "recvbuffer"); MPI_Alltoallv((void *)nodelist, scounts, sdispl, IDX_DATATYPE, (void *)recvbuffer, rcounts, rdispl, IDX_DATATYPE, *comm); /**************************************/ /* construct global node-element list */ /**************************************/ gnptr = idxsmalloc(my_nns+1, 0, "gnptr"); for (i=0; i<npes; i++) { for (j=rdispl[i]/2; j<rdispl[i+1]/2; j++) { lnode = recvbuffer[j].key-firstnode; ASSERTS(lnode >= 0 && lnode < my_nns) gnptr[lnode]++; } } MAKECSR(i, my_nns, gnptr); gnind = idxmalloc(amax(1, gnptr[my_nns]), "gnind"); for (pe=0; pe<npes; pe++) { firstelm = elmdist[pe]; for (j=rdispl[pe]/2; j<rdispl[pe+1]/2; j++) { lnode = recvbuffer[j].key-firstnode; gnind[gnptr[lnode]++] = recvbuffer[j].val+firstelm; } } SHIFTCSR(i, my_nns, gnptr); /*********************************************************/ /* send the node-element info to the relevant processors */ /*********************************************************/ iset(npes, 0, scounts); /* use a hash table to ensure that each node is sent to a proc only once */ for (pe=0; pe<npes; pe++) { for (j=rdispl[pe]/2; j<rdispl[pe+1]/2; j++) { lnode = recvbuffer[j].key-firstnode; if (htable[lnode] == -1) { scounts[pe] += gnptr[lnode+1]-gnptr[lnode]; htable[lnode] = 1; } } /* now reset the hash table */ for (j=rdispl[pe]/2; j<rdispl[pe+1]/2; j++) { lnode = recvbuffer[j].key-firstnode; htable[lnode] = -1; } } MPI_Alltoall((void *)scounts, 1, MPI_INT, (void *)rcounts, 1, MPI_INT, *comm); icopy(npes, scounts, sdispl); MAKECSR(i, npes, sdispl); /* create the send buffer */ nsend = sdispl[npes]; sbuffer = (idxtype *)realloc(nodelist, sizeof(idxtype)*amax(1, nsend)); count = 0; for (pe=0; pe<npes; pe++) { for (j=rdispl[pe]/2; j<rdispl[pe+1]/2; j++) { lnode = recvbuffer[j].key-firstnode; if (htable[lnode] == -1) { for (k=gnptr[lnode]; k<gnptr[lnode+1]; k++) { if (k == gnptr[lnode]) sbuffer[count++] = -1*(gnind[k]+1); else sbuffer[count++] = gnind[k]; } htable[lnode] = 1; } } ASSERTS(count == sdispl[pe+1]); /* now reset the hash table */ for (j=rdispl[pe]/2; j<rdispl[pe+1]/2; j++) { lnode = recvbuffer[j].key-firstnode; htable[lnode] = -1; } } icopy(npes, rcounts, rdispl); MAKECSR(i, npes, rdispl); nrecv = rdispl[npes]; rbuffer = (idxtype *)realloc(recvbuffer, sizeof(idxtype)*amax(1, nrecv)); MPI_Alltoallv((void *)sbuffer, scounts, sdispl, IDX_DATATYPE, (void *)rbuffer, rcounts, rdispl, IDX_DATATYPE, *comm); k = -1; nptr = idxsmalloc(lnns+1, 0, "nptr"); nind = rbuffer; for (pe=0; pe<npes; pe++) { for (j=rdispl[pe]; j<rdispl[pe+1]; j++) { if (nind[j] < 0) { k++; nind[j] = (-1*nind[j])-1; } nptr[k]++; } } MAKECSR(i, lnns, nptr); ASSERTS(k+1 == lnns); ASSERTS(nptr[lnns] == nrecv) myxadj = *xadj = idxsmalloc(nelms+1, 0, "xadj"); idxset(mask+1, -1, htable); firstelm = elmdist[mype]; /* Two passes -- in first pass, simply find out the memory requirements */ for (pass=0; pass<2; pass++) { for (i=0; i<nelms; i++) { for (count=0, j=eptr[i]; j<eptr[i+1]; j++) { node = eind[j]; for (k=nptr[node]; k<nptr[node+1]; k++) { if ((kk=nind[k]) == firstelm+i) continue; m = htable[(kk&mask)]; if (m == -1) { ind[count] = kk; wgt[count] = 1; htable[(kk&mask)] = count++; } else { if (ind[m] == kk) { wgt[m]++; } else { for (jj=0; jj<count; jj++) { if (ind[jj] == kk) { wgt[jj]++; break; } } if (jj == count) { ind[count] = kk; wgt[count++] = 1; } } } } } for (j=0; j<count; j++) { htable[(ind[j]&mask)] = -1; if (wgt[j] >= *ncommonnodes) { if (pass == 0) myxadj[i]++; else myadjncy[myxadj[i]++] = ind[j]; } } } if (pass == 0) { MAKECSR(i, nelms, myxadj); myadjncy = *adjncy = idxmalloc(myxadj[nelms], "adjncy"); } else { SHIFTCSR(i, nelms, myxadj); } } /*****************************************/ /* correctly renumber the elements array */ /*****************************************/ for (i=0; i<eptr[nelms]; i++) eind[i] = nmap[eind[i]] + gminnode; if (*numflag == 1) ChangeNumberingMesh2(elmdist, eptr, eind, myxadj, myadjncy, NULL, npes, mype, 0); /* do not free nodelist, recvbuffer, rbuffer */ GKfree((void **)&scounts, (void **)&nodedist, (void **)&nmap, (void **)&sbuffer, (void **)&htable, (void **)&nptr, (void **)&nind, (void **)&gnptr, (void **)&gnind, (void **)&auxarray, LTERM); FreeCtrl(&ctrl); return; }
/************************************************************************* * This function finds a matching **************************************************************************/ void Moc_GlobalMatch_Balance(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace) { int h, i, ii, j, k; int nnbrs, nvtxs, ncon, cnvtxs, firstvtx, lastvtx, maxi, maxidx, nkept; int otherlastvtx, nrequests, nchanged, pass, nmatched, wside; idxtype *xadj, *ladjncy, *adjwgt, *vtxdist, *home, *myhome, *shome, *rhome; idxtype *match, *rmatch, *smatch; idxtype *peind, *sendptr, *recvptr; idxtype *perm, *iperm, *nperm, *changed; floattype *nvwgt, maxnvwgt; int *nreqs_pe; KeyValueType *match_requests, *match_granted, *pe_requests; maxnvwgt = 1.0/((floattype)(ctrl->nparts)*MAXNVWGT_FACTOR); graph->match_type = MATCH_GLOBAL; IFSET(ctrl->dbglvl, DBG_TIME, MPI_Barrier(ctrl->comm)); IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->MatchTmr)); nvtxs = graph->nvtxs; ncon = graph->ncon; xadj = graph->xadj; ladjncy = graph->adjncy; adjwgt = graph->adjwgt; home = graph->home; nvwgt = graph->nvwgt; vtxdist = graph->vtxdist; firstvtx = vtxdist[ctrl->mype]; lastvtx = vtxdist[ctrl->mype+1]; match = graph->match = idxsmalloc(nvtxs+graph->nrecv, UNMATCHED, "HEM_Match: match"); myhome = idxsmalloc(nvtxs+graph->nrecv, UNMATCHED, "HEM_Match: myhome"); /*------------------------------------------------------------ / Send/Receive the home information of interface vertices /------------------------------------------------------------*/ if (ctrl->partType == ADAPTIVE_PARTITION || ctrl->partType == REFINE_PARTITION) { idxcopy(nvtxs, home, myhome); shome = wspace->indices; rhome = myhome + nvtxs; CommInterfaceData(ctrl, graph, myhome, shome, rhome); } nnbrs = graph->nnbrs; peind = graph->peind; sendptr = graph->sendptr; recvptr = graph->recvptr; /* Use wspace->indices as the tmp space for matching info of the boundary * vertices that are sent and received */ rmatch = match + nvtxs; smatch = wspace->indices; changed = smatch+graph->nsend; /* Use wspace->indices as the tmp space for match requests of the boundary * vertices that are sent and received */ match_requests = wspace->pairs; match_granted = match_requests + graph->nsend; nreqs_pe = ismalloc(nnbrs, 0, "Match_HEM: nreqs_pe"); nkept = graph->gnvtxs/ctrl->npes - nvtxs; perm = (idxtype *)wspace->degrees; iperm = perm + nvtxs; FastRandomPermute(nvtxs, perm, 1); for (i=0; i<nvtxs; i++) iperm[perm[i]] = i; nperm = iperm + nvtxs; for (i=0; i<nnbrs; i++) nperm[i] = i; /************************************************************* * Go now and find a matching by doing multiple iterations *************************************************************/ /* First nullify the heavy vertices */ for (nchanged=i=0; i<nvtxs; i++) { for (h=0; h<ncon; h++) if (nvwgt[i*ncon+h] > maxnvwgt) { break; } if (h != ncon) { match[i] = TOO_HEAVY; nchanged++; } } if (GlobalSESum(ctrl, nchanged) > 0) { IFSET(ctrl->dbglvl, DBG_PROGRESS, rprintf(ctrl, "We found %d heavy vertices!\n", GlobalSESum(ctrl, nchanged))); CommInterfaceData(ctrl, graph, match, smatch, rmatch); } for (nmatched=pass=0; pass<NMATCH_PASSES; pass++) { wside = (graph->level+pass)%2; nchanged = nrequests = 0; for (ii=nmatched; ii<nvtxs; ii++) { i = perm[ii]; if (match[i] == UNMATCHED) { /* Unmatched */ maxidx = i; maxi = -1; /* Find a heavy-edge matching */ for (j=xadj[i]; j<xadj[i+1]; j++) { k = ladjncy[j]; if (match[k] == UNMATCHED && myhome[k] == myhome[i] && (maxi == -1 || adjwgt[maxi] < adjwgt[j] || (maxidx < nvtxs && k < nvtxs && adjwgt[maxi] == adjwgt[j] && BetterVBalance(ncon,nvwgt+i*ncon,nvwgt+maxidx*ncon,nvwgt+k*ncon) >= 0))) { maxi = j; maxidx = k; } } if (maxi != -1) { k = ladjncy[maxi]; if (k < nvtxs) { /* Take care the local vertices first */ /* Here we give preference the local matching by granting it right away */ if (i <= k) { match[i] = firstvtx+k + KEEP_BIT; match[k] = firstvtx+i; } else { match[i] = firstvtx+k; match[k] = firstvtx+i + KEEP_BIT; } changed[nchanged++] = i; changed[nchanged++] = k; } else { /* Take care any remote boundary vertices */ match[k] = MAYBE_MATCHED; /* Alternate among which vertices will issue the requests */ if ((wside ==0 && firstvtx+i < graph->imap[k]) || (wside == 1 && firstvtx+i > graph->imap[k])) { match[i] = MAYBE_MATCHED; match_requests[nrequests].key = graph->imap[k]; match_requests[nrequests].val = firstvtx+i; nrequests++; } } } } } #ifdef DEBUG_MATCH PrintVector2(ctrl, nvtxs, firstvtx, match, "Match1"); myprintf(ctrl, "[c: %2d] Nlocal: %d, Nrequests: %d\n", c, nlocal, nrequests); #endif /*********************************************************** * Exchange the match_requests, requests for me are stored in * match_granted ************************************************************/ /* Issue the receives first. Note that from each PE can receive a maximum of the interface node that it needs to send it in the case of a mat-vec */ for (i=0; i<nnbrs; i++) { MPI_Irecv((void *)(match_granted+recvptr[i]), 2*(recvptr[i+1]-recvptr[i]), IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->rreq+i); } /* Issue the sends next. This needs some work */ ikeysort(nrequests, match_requests); for (j=i=0; i<nnbrs; i++) { otherlastvtx = vtxdist[peind[i]+1]; for (k=j; k<nrequests && match_requests[k].key < otherlastvtx; k++); MPI_Isend((void *)(match_requests+j), 2*(k-j), IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->sreq+i); j = k; } /* OK, now get into the loop waiting for the operations to finish */ MPI_Waitall(nnbrs, ctrl->rreq, ctrl->statuses); for (i=0; i<nnbrs; i++) { MPI_Get_count(ctrl->statuses+i, IDX_DATATYPE, nreqs_pe+i); nreqs_pe[i] = nreqs_pe[i]/2; /* Adjust for pairs of IDX_DATATYPE */ } MPI_Waitall(nnbrs, ctrl->sreq, ctrl->statuses); /*********************************************************** * Now, go and service the requests that you received in * match_granted ************************************************************/ RandomPermute(nnbrs, nperm, 0); for (ii=0; ii<nnbrs; ii++) { i = nperm[ii]; pe_requests = match_granted+recvptr[i]; for (j=0; j<nreqs_pe[i]; j++) { k = pe_requests[j].key; ASSERTP(ctrl, k >= firstvtx && k < lastvtx, (ctrl, "%d %d %d %d %d\n", firstvtx, lastvtx, k, j, peind[i])); /* myprintf(ctrl, "Requesting a match %d %d\n", pe_requests[j].key, pe_requests[j].val); */ if (match[k-firstvtx] == UNMATCHED) { /* Bingo, lets grant this request */ changed[nchanged++] = k-firstvtx; if (nkept >= 0) { /* Flip a coin for who gets it */ match[k-firstvtx] = pe_requests[j].val + KEEP_BIT; nkept--; } else { match[k-firstvtx] = pe_requests[j].val; pe_requests[j].key += KEEP_BIT; nkept++; } /* myprintf(ctrl, "Request from pe:%d (%d %d) granted!\n", peind[i], pe_requests[j].val, pe_requests[j].key); */ } else { /* We are not granting the request */ /* myprintf(ctrl, "Request from pe:%d (%d %d) not granted!\n", peind[i], pe_requests[j].val, pe_requests[j].key); */ pe_requests[j].key = UNMATCHED; } } } /*********************************************************** * Exchange the match_granted information. It is stored in * match_requests ************************************************************/ /* Issue the receives first. Note that from each PE can receive a maximum of the interface node that it needs to send during the case of a mat-vec */ for (i=0; i<nnbrs; i++) { MPI_Irecv((void *)(match_requests+sendptr[i]), 2*(sendptr[i+1]-sendptr[i]), IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->rreq+i); } /* Issue the sends next. */ for (i=0; i<nnbrs; i++) { MPI_Isend((void *)(match_granted+recvptr[i]), 2*nreqs_pe[i], IDX_DATATYPE, peind[i], 1, ctrl->comm, ctrl->sreq+i); } /* OK, now get into the loop waiting for the operations to finish */ MPI_Waitall(nnbrs, ctrl->rreq, ctrl->statuses); for (i=0; i<nnbrs; i++) { MPI_Get_count(ctrl->statuses+i, IDX_DATATYPE, nreqs_pe+i); nreqs_pe[i] = nreqs_pe[i]/2; /* Adjust for pairs of IDX_DATATYPE */ } MPI_Waitall(nnbrs, ctrl->sreq, ctrl->statuses); /*********************************************************** * Now, go and through the match_requests and update local * match information for the matchings that were granted. ************************************************************/ for (i=0; i<nnbrs; i++) { pe_requests = match_requests+sendptr[i]; for (j=0; j<nreqs_pe[i]; j++) { match[pe_requests[j].val-firstvtx] = pe_requests[j].key; if (pe_requests[j].key != UNMATCHED) changed[nchanged++] = pe_requests[j].val-firstvtx; } } for (i=0; i<nchanged; i++) { ii = iperm[changed[i]]; perm[ii] = perm[nmatched]; iperm[perm[nmatched]] = ii; nmatched++; } CommChangedInterfaceData(ctrl, graph, nchanged, changed, match, match_requests, match_granted, wspace->pv4); } /* Traverse the vertices and those that were unmatched, match them with themselves */ cnvtxs = 0; for (i=0; i<nvtxs; i++) { if (match[i] == UNMATCHED || match[i] == TOO_HEAVY) { match[i] = (firstvtx+i) + KEEP_BIT; cnvtxs++; } else if (match[i] >= KEEP_BIT) { /* A matched vertex which I get to keep */ cnvtxs++; } } if (ctrl->dbglvl&DBG_MATCHINFO) { PrintVector2(ctrl, nvtxs, firstvtx, match, "Match"); myprintf(ctrl, "Cnvtxs: %d\n", cnvtxs); rprintf(ctrl, "Done with matching...\n"); } GKfree((void **)(&myhome), (void **)(&nreqs_pe), LTERM); IFSET(ctrl->dbglvl, DBG_TIME, MPI_Barrier(ctrl->comm)); IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->MatchTmr)); IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->ContractTmr)); Moc_Global_CreateCoarseGraph(ctrl, graph, wspace, cnvtxs); IFSET(ctrl->dbglvl, DBG_TIME, MPI_Barrier(ctrl->comm)); IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->ContractTmr)); }
/************************************************************************* * This function sorts a distributed list of KeyValueType in increasing * order, and uses it to compute a partition. It uses samplesort. **************************************************************************/ void PartSort(CtrlType *ctrl, GraphType *graph, KeyValueType *elmnts, WorkSpaceType *wspace) { int i, j, k, nvtxs, nrecv, npes=ctrl->npes, mype=ctrl->mype, firstvtx, lastvtx; idxtype *scounts, *rcounts, *vtxdist, *perm; KeyValueType *relmnts, *mypicks, *allpicks; nvtxs = graph->nvtxs; vtxdist = graph->vtxdist; scounts = wspace->pv1; rcounts = wspace->pv2; /* Allocate memory for the splitters */ mypicks = (KeyValueType *)GKmalloc(sizeof(KeyValueType)*(npes+1), "ParSort: mypicks"); allpicks = (KeyValueType *)GKmalloc(sizeof(KeyValueType)*npes*npes, "ParSort: allpicks"); /* Sort the local elements */ ikeysort(nvtxs, elmnts); /* Select the local npes-1 equally spaced elements */ for (i=1; i<npes; i++) { mypicks[i-1].key = elmnts[i*(nvtxs/npes)].key; mypicks[i-1].val = elmnts[i*(nvtxs/npes)].val; } /* PrintPairs(ctrl, npes-1, mypicks, "Mypicks"); */ /* Gather the picks to all the processors */ MPI_Allgather((void *)mypicks, 2*(npes-1), IDX_DATATYPE, (void *)allpicks, 2*(npes-1), IDX_DATATYPE, ctrl->comm); /* PrintPairs(ctrl, npes*(npes-1), allpicks, "Allpicks"); */ /* Sort all the picks */ ikeyvalsort(npes*(npes-1), allpicks); /* PrintPairs(ctrl, npes*(npes-1), allpicks, "Allpicks"); */ /* Select the final splitters. Set the boundaries to simplify coding */ for (i=1; i<npes; i++) mypicks[i] = allpicks[i*(npes-1)]; mypicks[0].key = MIN_INT; mypicks[npes].key = MAX_INT; /* PrintPairs(ctrl, npes+1, mypicks, "Mypicks"); */ /* Compute the number of elements that belong to each bucket */ idxset(npes, 0, scounts); for (j=i=0; i<nvtxs; i++) { if (elmnts[i].key < mypicks[j+1].key || (elmnts[i].key == mypicks[j+1].key && elmnts[i].val < mypicks[j+1].val)) scounts[j]++; else scounts[++j]++; } MPI_Alltoall(scounts, 1, IDX_DATATYPE, rcounts, 1, IDX_DATATYPE, ctrl->comm); /* PrintVector(ctrl, npes, 0, scounts, "Scounts"); PrintVector(ctrl, npes, 0, rcounts, "Rcounts"); */ /* Allocate memory for sorted elements and receive them */ MAKECSR(i, npes, scounts); MAKECSR(i, npes, rcounts); nrecv = rcounts[npes]; if (wspace->nlarge >= nrecv) relmnts = (KeyValueType *)wspace->pairs; else relmnts = (KeyValueType *)GKmalloc(sizeof(KeyValueType)*nrecv, "ParSort: relmnts"); /* Issue the receives first */ for (i=0; i<npes; i++) MPI_Irecv((void *)(relmnts+rcounts[i]), 2*(rcounts[i+1]-rcounts[i]), IDX_DATATYPE, i, 1, ctrl->comm, ctrl->rreq+i); /* Issue the sends next */ for (i=0; i<npes; i++) MPI_Isend((void *)(elmnts+scounts[i]), 2*(scounts[i+1]-scounts[i]), IDX_DATATYPE, i, 1, ctrl->comm, ctrl->sreq+i); MPI_Waitall(npes, ctrl->rreq, ctrl->statuses); MPI_Waitall(npes, ctrl->sreq, ctrl->statuses); /* OK, now do the local sort of the relmnts. Use perm to keep track original order */ perm = idxmalloc(nrecv, "ParSort: perm"); for (i=0; i<nrecv; i++) { perm[i] = relmnts[i].val; relmnts[i].val = i; } ikeysort(nrecv, relmnts); /* Compute what needs to be shifted */ MPI_Scan((void *)(&nrecv), (void *)(&lastvtx), 1, MPI_INT, MPI_SUM, ctrl->comm); firstvtx = lastvtx-nrecv; /*myprintf(ctrl, "first, last: %d %d\n", firstvtx, lastvtx); */ for (j=0, i=0; i<npes; i++) { if (vtxdist[i+1] > firstvtx) { /* Found the first PE that is passed me */ if (vtxdist[i+1] >= lastvtx) { /* myprintf(ctrl, "Shifting %d elements to processor %d\n", lastvtx-firstvtx, i); */ for (k=0; k<lastvtx-firstvtx; k++, j++) relmnts[relmnts[j].val].key = i; } else { /* myprintf(ctrl, "Shifting %d elements to processor %d\n", vtxdist[i+1]-firstvtx, i); */ for (k=0; k<vtxdist[i+1]-firstvtx; k++, j++) relmnts[relmnts[j].val].key = i; firstvtx = vtxdist[i+1]; } } if (vtxdist[i+1] >= lastvtx) break; } /* Reverse the ordering on the relmnts[].val */ for (i=0; i<nrecv; i++) { ASSERTP(ctrl, relmnts[i].key>=0 && relmnts[i].key<npes, (ctrl, "%d %d\n", i, relmnts[i].key)); relmnts[i].val = perm[i]; } /* OK, now sent it back */ /* Issue the receives first */ for (i=0; i<npes; i++) MPI_Irecv((void *)(elmnts+scounts[i]), 2*(scounts[i+1]-scounts[i]), IDX_DATATYPE, i, 1, ctrl->comm, ctrl->rreq+i); /* Issue the sends next */ for (i=0; i<npes; i++) MPI_Isend((void *)(relmnts+rcounts[i]), 2*(rcounts[i+1]-rcounts[i]), IDX_DATATYPE, i, 1, ctrl->comm, ctrl->sreq+i); MPI_Waitall(npes, ctrl->rreq, ctrl->statuses); MPI_Waitall(npes, ctrl->sreq, ctrl->statuses); /* Construct a partition for the graph */ graph->where = idxmalloc(graph->nvtxs+graph->nrecv, "PartSort: graph->where"); firstvtx = vtxdist[mype]; for (i=0; i<nvtxs; i++) { ASSERTP(ctrl, elmnts[i].key>=0 && elmnts[i].key<npes, (ctrl, "%d %d\n", i, elmnts[i].key)); ASSERTP(ctrl, elmnts[i].val>=vtxdist[mype] && elmnts[i].val<vtxdist[mype+1], (ctrl, "%d %d %d %d\n", i, vtxdist[mype], vtxdist[mype+1], elmnts[i].val)); graph->where[elmnts[i].val-firstvtx] = elmnts[i].key; } GKfree((void **)&mypicks, (void **)&allpicks, (void **)&perm, LTERM); if (wspace->nlarge < nrecv) free(relmnts); }