00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #include <parmetislib.h>
00016
00017
00018
00019
00020
00021
00022
00023 void ParMETIS_V3_PartKway(idxtype *vtxdist, idxtype *xadj, idxtype *adjncy, idxtype *vwgt,
00024 idxtype *adjwgt, int *wgtflag, int *numflag, int *ncon, int *nparts,
00025 floattype *tpwgts, floattype *ubvec, int *options, int *edgecut, idxtype *part,
00026 MPI_Comm *comm)
00027 {
00028 int h, i;
00029 int nvtxs = -1, npes, mype;
00030 CtrlType ctrl;
00031 WorkSpaceType wspace;
00032 GraphType *graph;
00033 floattype avg, maximb, *mytpwgts;
00034 int moptions[10];
00035 int seed, dbglvl = 0;
00036 int iwgtflag, inumflag, incon, inparts, ioptions[10];
00037 floattype *itpwgts, iubvec[MAXNCON];
00038
00039 MPI_Comm_size(*comm, &npes);
00040 MPI_Comm_rank(*comm, &mype);
00041
00042
00043
00044
00045
00046 if (options != NULL && options[0] == 1)
00047 dbglvl = options[PMV3_OPTION_DBGLVL];
00048
00049 CheckInputs(STATIC_PARTITION, npes, dbglvl, wgtflag, &iwgtflag, numflag, &inumflag, ncon,
00050 &incon, nparts, &inparts, tpwgts, &itpwgts, ubvec, iubvec, NULL, NULL,
00051 options, ioptions, part, comm);
00052
00053
00054
00055
00056
00057 if (inparts <= 1) {
00058 idxset(vtxdist[mype+1]-vtxdist[mype], 0, part);
00059 *edgecut = 0;
00060 return;
00061 }
00062
00063
00064
00065
00066 if (npes == 1 && inparts > 1) {
00067 moptions[0] = 0;
00068 nvtxs = vtxdist[1];
00069
00070 if (incon == 1) {
00071 METIS_WPartGraphKway(&nvtxs, xadj, adjncy, vwgt, adjwgt, &iwgtflag, &inumflag,
00072 &inparts, itpwgts, moptions, edgecut, part);
00073 }
00074 else {
00075
00076 mytpwgts = fmalloc(inparts, "mytpwgts");
00077 for (i=0; i<inparts; i++)
00078 mytpwgts[i] = itpwgts[i*incon];
00079
00080 moptions[7] = -1;
00081 METIS_mCPartGraphRecursive2(&nvtxs, &incon, xadj, adjncy, vwgt, adjwgt, &iwgtflag,
00082 &inumflag, &inparts, mytpwgts, moptions, edgecut, part);
00083
00084 free(mytpwgts);
00085 }
00086
00087 return;
00088 }
00089
00090
00091 if (inumflag == 1)
00092 ChangeNumbering(vtxdist, xadj, adjncy, part, npes, mype, 1);
00093
00094
00095
00096
00097 if (ioptions[0] == 1) {
00098 dbglvl = ioptions[PMV3_OPTION_DBGLVL];
00099 seed = ioptions[PMV3_OPTION_SEED];
00100 }
00101 else {
00102 dbglvl = GLOBAL_DBGLVL;
00103 seed = GLOBAL_SEED;
00104 }
00105 SetUpCtrl(&ctrl, inparts, dbglvl, *comm);
00106 ctrl.CoarsenTo = amin(vtxdist[npes]+1, 25*incon*amax(npes, inparts));
00107 ctrl.seed = (seed == 0) ? mype : seed*mype;
00108 ctrl.sync = GlobalSEMax(&ctrl, seed);
00109 ctrl.partType = STATIC_PARTITION;
00110 ctrl.ps_relation = -1;
00111 ctrl.tpwgts = itpwgts;
00112 scopy(incon, iubvec, ctrl.ubvec);
00113
00114 graph = Moc_SetUpGraph(&ctrl, incon, vtxdist, xadj, vwgt, adjncy, adjwgt, &iwgtflag);
00115
00116 PreAllocateMemory(&ctrl, graph, &wspace);
00117
00118 IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl));
00119 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00120 IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr));
00121
00122
00123
00124
00125
00126
00127
00128
00129 if (vtxdist[npes] < SMALLGRAPH || vtxdist[npes] < npes*20 || GlobalSESum(&ctrl, graph->nedges) == 0) {
00130 IFSET(ctrl.dbglvl, DBG_INFO, rprintf(&ctrl, "Partitioning a graph of size %d serially\n", vtxdist[npes]));
00131 PartitionSmallGraph(&ctrl, graph, &wspace);
00132 }
00133 else {
00134
00135
00136
00137 Moc_Global_Partition(&ctrl, graph, &wspace);
00138 ParallelReMapGraph(&ctrl, graph, &wspace);
00139 }
00140
00141 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00142 IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr));
00143
00144 idxcopy(graph->nvtxs, graph->where, part);
00145 *edgecut = graph->mincut;
00146
00147
00148
00149
00150 IFSET(ctrl.dbglvl, DBG_TIME, PrintTimingInfo(&ctrl));
00151 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00152
00153 if (ctrl.dbglvl&DBG_INFO) {
00154 rprintf(&ctrl, "Final %d-way CUT: %6d \tBalance: ", inparts, graph->mincut);
00155 avg = 0.0;
00156 for (h=0; h<incon; h++) {
00157 maximb = 0.0;
00158 for (i=0; i<inparts; i++)
00159 maximb = amax(maximb, graph->gnpwgts[i*incon+h]/itpwgts[i*incon+h]);
00160 avg += maximb;
00161 rprintf(&ctrl, "%.3f ", maximb);
00162 }
00163 rprintf(&ctrl, " avg: %.3f\n", avg/(floattype)incon);
00164 }
00165
00166 GKfree((void **)&itpwgts, (void **)&graph->lnpwgts, (void **)&graph->gnpwgts, (void **)&graph->nvwgt, LTERM);
00167 FreeInitialGraphAndRemap(graph, iwgtflag);
00168 FreeWSpace(&wspace);
00169 FreeCtrl(&ctrl);
00170
00171 if (inumflag == 1)
00172 ChangeNumbering(vtxdist, xadj, adjncy, part, npes, mype, 0);
00173
00174 }
00175
00176
00177
00178
00179
00180
00181 void Moc_Global_Partition(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00182 {
00183 int i, ncon, nparts;
00184 floattype ftmp, ubavg, lbavg, lbvec[MAXNCON];
00185
00186 ncon = graph->ncon;
00187 nparts = ctrl->nparts;
00188 ubavg = savg(graph->ncon, ctrl->ubvec);
00189
00190 SetUp(ctrl, graph, wspace);
00191
00192 if (ctrl->dbglvl&DBG_PROGRESS) {
00193 rprintf(ctrl, "[%6d %8d %5d %5d] [%d] [", graph->gnvtxs, GlobalSESum(ctrl, graph->nedges),
00194 GlobalSEMin(ctrl, graph->nvtxs), GlobalSEMax(ctrl, graph->nvtxs), ctrl->CoarsenTo);
00195 for (i=0; i<ncon; i++)
00196 rprintf(ctrl, " %.3f", GlobalSEMinFloat(ctrl,graph->nvwgt[samin_strd(graph->nvtxs, graph->nvwgt+i, ncon)*ncon+i]));
00197 rprintf(ctrl, "] [");
00198 for (i=0; i<ncon; i++)
00199 rprintf(ctrl, " %.3f", GlobalSEMaxFloat(ctrl, graph->nvwgt[samax_strd(graph->nvtxs, graph->nvwgt+i, ncon)*ncon+i]));
00200 rprintf(ctrl, "]\n");
00201 }
00202
00203 if (graph->gnvtxs < 1.3*ctrl->CoarsenTo ||
00204 (graph->finer != NULL &&
00205 graph->gnvtxs > graph->finer->gnvtxs*COARSEN_FRACTION)) {
00206
00207
00208 graph->where = idxmalloc(graph->nvtxs+graph->nrecv, "graph->where");
00209 Moc_InitPartition_RB(ctrl, graph, wspace);
00210
00211 if (ctrl->dbglvl&DBG_PROGRESS) {
00212 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00213 rprintf(ctrl, "nvtxs: %10d, balance: ", graph->gnvtxs);
00214 for (i=0; i<graph->ncon; i++)
00215 rprintf(ctrl, "%.3f ", lbvec[i]);
00216 rprintf(ctrl, "\n");
00217 }
00218
00219
00220 if (graph->finer == NULL) {
00221 Moc_ComputePartitionParams(ctrl, graph, wspace);
00222 Moc_KWayFM(ctrl, graph, wspace, NGR_PASSES);
00223 }
00224 }
00225 else {
00226 Moc_GlobalMatch_Balance(ctrl, graph, wspace);
00227
00228 Moc_Global_Partition(ctrl, graph->coarser, wspace);
00229
00230 Moc_ProjectPartition(ctrl, graph, wspace);
00231 Moc_ComputePartitionParams(ctrl, graph, wspace);
00232
00233 if (graph->ncon > 1 && graph->level < 3) {
00234 for (i=0; i<ncon; i++) {
00235 ftmp = ssum_strd(nparts, graph->gnpwgts+i, ncon);
00236 if (ftmp != 0.0)
00237 lbvec[i] = (floattype)(nparts) *
00238 graph->gnpwgts[samax_strd(nparts, graph->gnpwgts+i, ncon)*ncon+i]/ftmp;
00239 else
00240 lbvec[i] = 1.0;
00241 }
00242 lbavg = savg(graph->ncon, lbvec);
00243
00244 if (lbavg > ubavg + 0.035) {
00245 if (ctrl->dbglvl&DBG_PROGRESS) {
00246 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00247 rprintf(ctrl, "nvtxs: %10d, cut: %8d, balance: ", graph->gnvtxs, graph->mincut);
00248 for (i=0; i<graph->ncon; i++)
00249 rprintf(ctrl, "%.3f ", lbvec[i]);
00250 rprintf(ctrl, "\n");
00251 }
00252
00253 Moc_KWayBalance(ctrl, graph, wspace, graph->ncon);
00254 }
00255 }
00256
00257 Moc_KWayFM(ctrl, graph, wspace, NGR_PASSES);
00258
00259 if (ctrl->dbglvl&DBG_PROGRESS) {
00260 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00261 rprintf(ctrl, "nvtxs: %10d, cut: %8d, balance: ", graph->gnvtxs, graph->mincut);
00262 for (i=0; i<graph->ncon; i++)
00263 rprintf(ctrl, "%.3f ", lbvec[i]);
00264 rprintf(ctrl, "\n");
00265 }
00266
00267 if (graph->level != 0)
00268 GKfree((void **)&graph->lnpwgts, (void **)&graph->gnpwgts, LTERM);
00269 }
00270
00271 return;
00272 }
00273
00274