00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include <parmetislib.h>
00015
00016
00017
00018
00019
00020
00021
00022 void Balance_Partition(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00023 {
00024 int i, j, mype, npes, nvtxs, nedges, ncon;
00025 idxtype *vtxdist, *xadj, *adjncy, *adjwgt, *vwgt, *vsize;
00026 idxtype *part, *lwhere, *home;
00027 GraphType *agraph, cgraph;
00028 CtrlType myctrl;
00029 int lnparts, fpart, fpe, lnpes, ngroups, srnpes, srmype;
00030 int twoparts=2, numflag = 0, wgtflag = 3, moptions[10], edgecut, max_cut;
00031 int sr_pe, gd_pe, sr, gd, who_wins, *rcounts, *rdispls;
00032 floattype my_cut, my_totalv, my_cost = -1.0, my_balance = -1.0, wsum;
00033 floattype rating, max_rating, your_cost = -1.0, your_balance = -1.0;
00034 floattype lbvec[MAXNCON], lbsum, min_lbsum, *mytpwgts, mytpwgts2[2], buffer[2];
00035 MPI_Status status;
00036 MPI_Comm ipcomm, srcomm;
00037 struct {
00038 floattype cost;
00039 int rank;
00040 } lpecost, gpecost;
00041
00042 IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->InitPartTmr));
00043
00044 vtxdist = graph->vtxdist;
00045 agraph = Moc_AssembleAdaptiveGraph(ctrl, graph, wspace);
00046 nvtxs = cgraph.nvtxs = agraph->nvtxs;
00047 nedges = cgraph.nedges = agraph->nedges;
00048 ncon = cgraph.ncon = agraph->ncon;
00049
00050 xadj = cgraph.xadj = idxmalloc(nvtxs*(5+ncon)+1+nedges*2, "U_IP: xadj");
00051 vwgt = cgraph.vwgt = xadj + nvtxs+1;
00052 vsize = cgraph.vsize = xadj + nvtxs*(1+ncon)+1;
00053 cgraph.where = agraph->where = part = xadj + nvtxs*(2+ncon)+1;
00054 lwhere = xadj + nvtxs*(3+ncon)+1;
00055 home = xadj + nvtxs*(4+ncon)+1;
00056 adjncy = cgraph.adjncy = xadj + nvtxs*(5+ncon)+1;
00057 adjwgt = cgraph.adjwgt = xadj + nvtxs*(5+ncon)+1 + nedges;
00058
00059
00060
00061 mytpwgts = fsmalloc(ctrl->nparts, 0.0, "mytpwgts");
00062 for (i=0; i<ctrl->nparts; i++)
00063 for (j=0; j<ncon; j++)
00064 mytpwgts[i] += ctrl->tpwgts[i*ncon+j];
00065 for (i=0; i<ctrl->nparts; i++)
00066 mytpwgts[i] /= (floattype)ncon;
00067
00068 idxcopy(nvtxs+1, agraph->xadj, xadj);
00069 idxcopy(nvtxs*ncon, agraph->vwgt, vwgt);
00070 idxcopy(nvtxs, agraph->vsize, vsize);
00071 idxcopy(nedges, agraph->adjncy, adjncy);
00072 idxcopy(nedges, agraph->adjwgt, adjwgt);
00073
00074
00075
00076 if (ctrl->ps_relation == DISCOUPLED) {
00077 rcounts = imalloc(ctrl->npes, "rcounts");
00078 rdispls = imalloc(ctrl->npes+1, "rdispls");
00079
00080 for (i=0; i<ctrl->npes; i++) {
00081 rdispls[i] = rcounts[i] = vtxdist[i+1]-vtxdist[i];
00082 }
00083 MAKECSR(i, ctrl->npes, rdispls);
00084
00085 MPI_Allgatherv((void *)graph->home, graph->nvtxs, IDX_DATATYPE,
00086 (void *)part, rcounts, rdispls, IDX_DATATYPE, ctrl->comm);
00087
00088 for (i=0; i<agraph->nvtxs; i++)
00089 home[i] = part[i];
00090
00091 GKfree((void **)&rcounts, (void **)&rdispls, LTERM);
00092 }
00093 else {
00094 for (i=0; i<ctrl->npes; i++)
00095 for (j=vtxdist[i]; j<vtxdist[i+1]; j++)
00096 part[j] = home[j] = i;
00097 }
00098
00099
00100 for (i=0; i<agraph->nvtxs; i++) {
00101 if (part[i] >= ctrl->nparts)
00102 part[i] = home[i] = part[i] % ctrl->nparts;
00103 if (part[i] < 0)
00104 part[i] = home[i] = (-1*part[i]) % ctrl->nparts;
00105 }
00106
00107
00108
00109 IFSET(ctrl->dbglvl, DBG_REFINEINFO, Moc_ComputeSerialBalance(ctrl, agraph, agraph->where, lbvec));
00110 IFSET(ctrl->dbglvl, DBG_REFINEINFO, rprintf(ctrl, "input cut: %d, balance: ", ComputeSerialEdgeCut(agraph)));
00111 for (i=0; i<agraph->ncon; i++)
00112 IFSET(ctrl->dbglvl, DBG_REFINEINFO, rprintf(ctrl, "%.3f ", lbvec[i]));
00113 IFSET(ctrl->dbglvl, DBG_REFINEINFO, rprintf(ctrl, "\n"));
00114
00115
00116
00117
00118 sr = (ctrl->mype % 2 == 0) ? 1 : 0;
00119 gd = (ctrl->mype % 2 == 1) ? 1 : 0;
00120
00121 if (graph->ncon > MAX_NCON_FOR_DIFFUSION || ctrl->npes == 1) {
00122 sr = 1;
00123 gd = 0;
00124 }
00125
00126 sr_pe = 0;
00127 gd_pe = 1;
00128
00129 MPI_Comm_split(ctrl->gcomm, sr, 0, &ipcomm);
00130 MPI_Comm_rank(ipcomm, &mype);
00131 MPI_Comm_size(ipcomm, &npes);
00132
00133 myctrl.dbglvl = 0;
00134 myctrl.mype = mype;
00135 myctrl.npes = npes;
00136 myctrl.comm = ipcomm;
00137 myctrl.sync = ctrl->sync;
00138 myctrl.seed = ctrl->seed;
00139 myctrl.nparts = ctrl->nparts;
00140 myctrl.ipc_factor = ctrl->ipc_factor;
00141 myctrl.redist_factor = ctrl->redist_base;
00142 myctrl.partType = ADAPTIVE_PARTITION;
00143 myctrl.ps_relation = DISCOUPLED;
00144 myctrl.tpwgts = ctrl->tpwgts;
00145 icopy(ncon, ctrl->tvwgts, myctrl.tvwgts);
00146 icopy(ncon, ctrl->ubvec, myctrl.ubvec);
00147
00148 if (sr == 1) {
00149
00150
00151
00152 ngroups = amax(amin(RIP_SPLIT_FACTOR, npes), 1);
00153 MPI_Comm_split(ipcomm, mype % ngroups, 0, &srcomm);
00154 MPI_Comm_rank(srcomm, &srmype);
00155 MPI_Comm_size(srcomm, &srnpes);
00156
00157 moptions[0] = 0;
00158 moptions[7] = ctrl->sync + (mype % ngroups) + 1;
00159
00160 idxset(nvtxs, 0, lwhere);
00161 lnparts = ctrl->nparts;
00162 fpart = fpe = 0;
00163 lnpes = srnpes;
00164 while (lnpes > 1 && lnparts > 1) {
00165 ASSERT(ctrl, agraph->nvtxs > 1);
00166
00167 mytpwgts2[0] = ssum(lnparts/2, mytpwgts+fpart);
00168 mytpwgts2[1] = 1.0-mytpwgts2[0];
00169
00170
00171 if (agraph->ncon == 1) {
00172 METIS_WPartGraphKway2(&agraph->nvtxs, agraph->xadj, agraph->adjncy, agraph->vwgt,
00173 agraph->adjwgt, &wgtflag, &numflag, &twoparts, mytpwgts2, moptions, &edgecut,
00174 part);
00175 }
00176 else {
00177 METIS_mCPartGraphRecursive2(&agraph->nvtxs, &ncon, agraph->xadj, agraph->adjncy,
00178 agraph->vwgt, agraph->adjwgt, &wgtflag, &numflag, &twoparts, mytpwgts2,
00179 moptions, &edgecut, part);
00180 }
00181
00182 wsum = ssum(lnparts/2, mytpwgts+fpart);
00183 sscale(lnparts/2, 1.0/wsum, mytpwgts+fpart);
00184 sscale(lnparts-lnparts/2, 1.0/(1.0-wsum), mytpwgts+fpart+lnparts/2);
00185
00186
00187 if (srmype < fpe+lnpes/2) {
00188 Moc_KeepPart(agraph, wspace, part, 0);
00189 lnpes = lnpes/2;
00190 lnparts = lnparts/2;
00191 }
00192 else {
00193 Moc_KeepPart(agraph, wspace, part, 1);
00194 fpart = fpart + lnparts/2;
00195 fpe = fpe + lnpes/2;
00196 lnpes = lnpes - lnpes/2;
00197 lnparts = lnparts - lnparts/2;
00198 }
00199 }
00200
00201
00202 if (lnparts == 1) {
00203
00204 if (srmype == fpe) {
00205 for (i=0; i<agraph->nvtxs; i++)
00206 lwhere[agraph->label[i]] = fpart;
00207 }
00208 }
00209
00210 else {
00211 if (ncon == 1)
00212 METIS_WPartGraphKway2(&agraph->nvtxs, agraph->xadj, agraph->adjncy, agraph->vwgt,
00213 agraph->adjwgt, &wgtflag, &numflag, &lnparts, mytpwgts+fpart, moptions,
00214 &edgecut, part);
00215 else
00216 METIS_mCPartGraphRecursive2(&agraph->nvtxs, &ncon, agraph->xadj, agraph->adjncy,
00217 agraph->vwgt, agraph->adjwgt, &wgtflag, &numflag, &lnparts, mytpwgts+fpart,
00218 moptions, &edgecut, part);
00219
00220 for (i=0; i<agraph->nvtxs; i++)
00221 lwhere[agraph->label[i]] = fpart + part[i];
00222 }
00223
00224 MPI_Allreduce((void *)lwhere, (void *)part, nvtxs, IDX_DATATYPE, MPI_SUM, srcomm);
00225
00226 edgecut = ComputeSerialEdgeCut(&cgraph);
00227 Moc_ComputeSerialBalance(ctrl, &cgraph, part, lbvec);
00228 lbsum = ssum(ncon, lbvec);
00229 MPI_Allreduce((void *)&edgecut, (void *)&max_cut, 1, MPI_INT, MPI_MAX, ipcomm);
00230 MPI_Allreduce((void *)&lbsum, (void *)&min_lbsum, 1, MPI_DOUBLE, MPI_MIN, ipcomm);
00231 lpecost.rank = ctrl->mype;
00232 lpecost.cost = lbsum;
00233 if (min_lbsum < UNBALANCE_FRACTION * (floattype)(ncon)) {
00234 if (lbsum < UNBALANCE_FRACTION * (floattype)(ncon))
00235 lpecost.cost = (floattype)edgecut;
00236 else
00237 lpecost.cost = (floattype)max_cut + lbsum;
00238 }
00239 MPI_Allreduce((void *)&lpecost, (void *)&gpecost, 1, MPI_DOUBLE_INT, MPI_MINLOC, ipcomm);
00240
00241 if (ctrl->mype == gpecost.rank && ctrl->mype != sr_pe) {
00242 MPI_Send((void *)part, nvtxs, IDX_DATATYPE, sr_pe, 1, ctrl->comm);
00243 }
00244
00245 if (ctrl->mype != gpecost.rank && ctrl->mype == sr_pe) {
00246 MPI_Recv((void *)part, nvtxs, IDX_DATATYPE, gpecost.rank, 1, ctrl->comm, &status);
00247 }
00248
00249 if (ctrl->mype == sr_pe) {
00250 idxcopy(nvtxs, part, lwhere);
00251 SerialRemap(&cgraph, ctrl->nparts, home, lwhere, part, ctrl->tpwgts);
00252 }
00253
00254 MPI_Comm_free(&srcomm);
00255 }
00256
00257
00258
00259 else {
00260
00261
00262
00263 MPI_Comm_split(ipcomm, MPI_UNDEFINED, 0, &srcomm);
00264
00265 if (ncon == 1) {
00266 rating = WavefrontDiffusion(&myctrl, agraph, home);
00267 Moc_ComputeSerialBalance(ctrl, &cgraph, part, lbvec);
00268 lbsum = ssum(ncon, lbvec);
00269
00270
00271 MPI_Allreduce((void *)&rating, (void *)&max_rating, 1, MPI_DOUBLE, MPI_MAX, ipcomm);
00272 MPI_Allreduce((void *)&lbsum, (void *)&min_lbsum, 1, MPI_DOUBLE, MPI_MIN, ipcomm);
00273
00274 lpecost.rank = ctrl->mype;
00275 lpecost.cost = lbsum;
00276 if (min_lbsum < UNBALANCE_FRACTION * (floattype)(ncon)) {
00277 if (lbsum < UNBALANCE_FRACTION * (floattype)(ncon))
00278 lpecost.cost = rating;
00279 else
00280 lpecost.cost = max_rating + lbsum;
00281 }
00282
00283 MPI_Allreduce((void *)&lpecost, (void *)&gpecost, 1, MPI_DOUBLE_INT, MPI_MINLOC, ipcomm);
00284
00285
00286 if (ctrl->mype == gpecost.rank && ctrl->mype != gd_pe)
00287 MPI_Send((void *)part, nvtxs, IDX_DATATYPE, gd_pe, 1, ctrl->comm);
00288
00289 if (ctrl->mype != gpecost.rank && ctrl->mype == gd_pe)
00290 MPI_Recv((void *)part, nvtxs, IDX_DATATYPE, gpecost.rank, 1, ctrl->comm, &status);
00291
00292 if (ctrl->mype == gd_pe) {
00293 idxcopy(nvtxs, part, lwhere);
00294 SerialRemap(&cgraph, ctrl->nparts, home, lwhere, part, ctrl->tpwgts);
00295 }
00296 }
00297 else {
00298 Moc_Diffusion(&myctrl, agraph, graph->vtxdist, agraph->where, home, wspace, N_MOC_GD_PASSES);
00299 }
00300 }
00301
00302 if (graph->ncon <= MAX_NCON_FOR_DIFFUSION) {
00303 if (ctrl->mype == sr_pe || ctrl->mype == gd_pe) {
00304
00305
00306
00307 my_cut = (floattype) ComputeSerialEdgeCut(&cgraph);
00308 my_totalv = (floattype) Mc_ComputeSerialTotalV(&cgraph, home);
00309 Moc_ComputeSerialBalance(ctrl, &cgraph, part, lbvec);
00310 my_balance = ssum(cgraph.ncon, lbvec);
00311 my_balance /= (floattype) cgraph.ncon;
00312 my_cost = ctrl->ipc_factor * my_cut + REDIST_WGT * ctrl->redist_base * my_totalv;
00313
00314 IFSET(ctrl->dbglvl, DBG_REFINEINFO, printf("%s initial cut: %.1f, totalv: %.1f, balance: %.3f\n",
00315 (ctrl->mype == sr_pe ? "scratch-remap" : "diffusion"), my_cut, my_totalv, my_balance));
00316
00317 if (ctrl->mype == gd_pe) {
00318 buffer[0] = my_cost;
00319 buffer[1] = my_balance;
00320 MPI_Send((void *)buffer, 2, MPI_DOUBLE, sr_pe, 1, ctrl->comm);
00321 }
00322 else {
00323 MPI_Recv((void *)buffer, 2, MPI_DOUBLE, gd_pe, 1, ctrl->comm, &status);
00324 your_cost = buffer[0];
00325 your_balance = buffer[1];
00326 }
00327 }
00328
00329 if (ctrl->mype == sr_pe) {
00330 who_wins = gd_pe;
00331 if ((my_balance < 1.1 && your_balance > 1.1) ||
00332 (my_balance < 1.1 && your_balance < 1.1 && my_cost < your_cost) ||
00333 (my_balance > 1.1 && your_balance > 1.1 && my_balance < your_balance)) {
00334 who_wins = sr_pe;
00335 }
00336 }
00337
00338 MPI_Bcast((void *)&who_wins, 1, MPI_INT, sr_pe, ctrl->comm);
00339 }
00340 else {
00341 who_wins = sr_pe;
00342 }
00343
00344 MPI_Bcast((void *)part, nvtxs, IDX_DATATYPE, who_wins, ctrl->comm);
00345 idxcopy(graph->nvtxs, part+vtxdist[ctrl->mype], graph->where);
00346
00347 MPI_Comm_free(&ipcomm);
00348 GKfree((void **)&xadj, (void **)&mytpwgts, LTERM);
00349
00350
00351
00352
00353 GKfree((void **)&agraph->xadj, (void **)&agraph->adjncy, (void **)&agraph->vwgt, (void **)&agraph->nvwgt, LTERM);
00354 GKfree((void **)&agraph->vsize, (void **)&agraph->adjwgt, (void **)&agraph->label, LTERM);
00355 GKfree((void **)&agraph, LTERM);
00356
00357 IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->InitPartTmr));
00358
00359 }
00360
00361
00362
00363
00364
00365
00366 GraphType *Moc_AssembleAdaptiveGraph(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00367 {
00368 int i, j, k, l, gnvtxs, nvtxs, ncon, gnedges, nedges, gsize;
00369 idxtype *xadj, *vwgt, *vsize, *adjncy, *adjwgt, *vtxdist, *imap;
00370 idxtype *axadj, *aadjncy, *aadjwgt, *avwgt, *avsize = NULL, *alabel;
00371 idxtype *mygraph, *ggraph;
00372 int *rcounts, *rdispls, mysize;
00373 floattype *anvwgt;
00374 GraphType *agraph;
00375
00376 gnvtxs = graph->gnvtxs;
00377 nvtxs = graph->nvtxs;
00378 ncon = graph->ncon;
00379 nedges = graph->xadj[nvtxs];
00380 xadj = graph->xadj;
00381 vwgt = graph->vwgt;
00382 vsize = graph->vsize;
00383 adjncy = graph->adjncy;
00384 adjwgt = graph->adjwgt;
00385 vtxdist = graph->vtxdist;
00386 imap = graph->imap;
00387
00388
00389
00390
00391 rcounts = imalloc(ctrl->npes, "AssembleGraph: rcounts");
00392 switch (ctrl->partType) {
00393 case STATIC_PARTITION:
00394 mysize = (1+ncon)*nvtxs + 2*nedges;
00395 break;
00396 case ADAPTIVE_PARTITION:
00397 case REFINE_PARTITION:
00398 mysize = (2+ncon)*nvtxs + 2*nedges;
00399 break;
00400 default:
00401 printf("WARNING: bad value for ctrl->partType %d\n", ctrl->partType);
00402 break;
00403 }
00404 MPI_Allgather((void *)(&mysize), 1, MPI_INT, (void *)rcounts, 1, MPI_INT, ctrl->comm);
00405
00406 rdispls = imalloc(ctrl->npes+1, "AssembleGraph: rdispls");
00407 rdispls[0] = 0;
00408 for (i=1; i<ctrl->npes+1; i++)
00409 rdispls[i] = rdispls[i-1] + rcounts[i-1];
00410
00411
00412 mygraph = (mysize <= wspace->maxcore ? wspace->core : idxmalloc(mysize, "AssembleGraph: mygraph"));
00413 for (k=i=0; i<nvtxs; i++) {
00414 mygraph[k++] = xadj[i+1]-xadj[i];
00415 for (j=0; j<ncon; j++)
00416 mygraph[k++] = vwgt[i*ncon+j];
00417 if (ctrl->partType == ADAPTIVE_PARTITION || ctrl->partType == REFINE_PARTITION)
00418 mygraph[k++] = vsize[i];
00419 for (j=xadj[i]; j<xadj[i+1]; j++) {
00420 mygraph[k++] = imap[adjncy[j]];
00421 mygraph[k++] = adjwgt[j];
00422 }
00423 }
00424 ASSERT(ctrl, mysize == k);
00425
00426
00427
00428
00429 gsize = rdispls[ctrl->npes];
00430 ggraph = (gsize <= wspace->maxcore-mysize ? wspace->core+mysize : idxmalloc(gsize, "AssembleGraph: ggraph"));
00431 MPI_Allgatherv((void *)mygraph, mysize, IDX_DATATYPE, (void *)ggraph, rcounts, rdispls, IDX_DATATYPE, ctrl->comm);
00432
00433 GKfree((void **)&rcounts, (void **)&rdispls, LTERM);
00434 if (mysize > wspace->maxcore)
00435 free(mygraph);
00436
00437 agraph = CreateGraph();
00438 agraph->nvtxs = gnvtxs;
00439 switch (ctrl->partType) {
00440 case STATIC_PARTITION:
00441 agraph->nedges = gnedges = (gsize-(1+ncon)*gnvtxs)/2;
00442 break;
00443 case ADAPTIVE_PARTITION:
00444 case REFINE_PARTITION:
00445 agraph->nedges = gnedges = (gsize-(2+ncon)*gnvtxs)/2;
00446 break;
00447 default:
00448 printf("WARNING: bad value for ctrl->partType %d\n", ctrl->partType);
00449 agraph->nedges = gnedges = -1;
00450 break;
00451 }
00452
00453 agraph->ncon = ncon;
00454
00455
00456
00457
00458 axadj = agraph->xadj = idxmalloc(gnvtxs+1, "AssembleGraph: axadj");
00459 avwgt = agraph->vwgt = idxmalloc(gnvtxs*ncon, "AssembleGraph: avwgt");
00460 anvwgt = agraph->nvwgt = fmalloc(gnvtxs*ncon, "AssembleGraph: anvwgt");
00461 aadjncy = agraph->adjncy = idxmalloc(gnedges, "AssembleGraph: adjncy");
00462 aadjwgt = agraph->adjwgt = idxmalloc(gnedges, "AssembleGraph: adjwgt");
00463 alabel = agraph->label = idxmalloc(gnvtxs, "AssembleGraph: alabel");
00464 if (ctrl->partType == ADAPTIVE_PARTITION || ctrl->partType == REFINE_PARTITION)
00465 avsize = agraph->vsize = idxmalloc(gnvtxs, "AssembleGraph: avsize");
00466
00467 for (k=j=i=0; i<gnvtxs; i++) {
00468 axadj[i] = ggraph[k++];
00469 for (l=0; l<ncon; l++)
00470 avwgt[i*ncon+l] = ggraph[k++];
00471 if (ctrl->partType == ADAPTIVE_PARTITION || ctrl->partType == REFINE_PARTITION)
00472 avsize[i] = ggraph[k++];
00473 for (l=0; l<axadj[i]; l++) {
00474 aadjncy[j] = ggraph[k++];
00475 aadjwgt[j] = ggraph[k++];
00476 j++;
00477 }
00478 }
00479
00480
00481
00482
00483 MAKECSR(i, gnvtxs, axadj);
00484
00485 for (i=0; i<gnvtxs; i++)
00486 for (j=0; j<ncon; j++)
00487 anvwgt[i*ncon+j] = (floattype)(agraph->vwgt[i*ncon+j]) / (floattype)(ctrl->tvwgts[j]);
00488
00489 for (i=0; i<gnvtxs; i++)
00490 alabel[i] = i;
00491
00492 if (gsize > wspace->maxcore-mysize)
00493 free(ggraph);
00494
00495 return agraph;
00496 }
00497
00498