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
00024 void ParMETIS_V3_AdaptiveRepart(idxtype *vtxdist, idxtype *xadj, idxtype *adjncy,
00025 idxtype *vwgt, idxtype *vsize, idxtype *adjwgt, int *wgtflag, int *numflag,
00026 int *ncon, int *nparts, floattype *tpwgts, floattype *ubvec, floattype *ipc2redist,
00027 int *options, int *edgecut, idxtype *part, MPI_Comm *comm)
00028 {
00029 int h, i;
00030 int npes, mype;
00031 CtrlType ctrl;
00032 WorkSpaceType wspace;
00033 GraphType *graph;
00034 int tewgt, tvsize, nmoved, maxin, maxout, vtx_factor;
00035 floattype gtewgt, gtvsize, avg, maximb;
00036 int ps_relation, seed, dbglvl = 0;
00037 int iwgtflag, inumflag, incon, inparts, ioptions[10];
00038 floattype iipc2redist, *itpwgts, iubvec[MAXNCON];
00039
00040 MPI_Comm_size(*comm, &npes);
00041 MPI_Comm_rank(*comm, &mype);
00042
00043
00044
00045
00046 if (options != NULL && options[0] == 1)
00047 dbglvl = options[PMV3_OPTION_DBGLVL];
00048 CheckInputs(ADAPTIVE_PARTITION, npes, dbglvl, wgtflag, &iwgtflag, numflag, &inumflag,
00049 ncon, &incon, nparts, &inparts, tpwgts, &itpwgts, ubvec, iubvec,
00050 ipc2redist, &iipc2redist, options, ioptions, part, comm);
00051
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 (inumflag == 1)
00067 ChangeNumbering(vtxdist, xadj, adjncy, part, npes, mype, 1);
00068
00069
00070
00071
00072 if (ioptions[0] == 1) {
00073 dbglvl = ioptions[PMV3_OPTION_DBGLVL];
00074 seed = ioptions[PMV3_OPTION_SEED];
00075 ps_relation = (npes == inparts ? ioptions[PMV3_OPTION_PSR] : DISCOUPLED);
00076 }
00077 else {
00078 dbglvl = GLOBAL_DBGLVL;
00079 seed = GLOBAL_SEED;
00080 ps_relation = (npes == inparts ? COUPLED : DISCOUPLED);
00081 }
00082
00083 SetUpCtrl(&ctrl, inparts, dbglvl, *comm);
00084 vtx_factor = (amax(npes, inparts) > 256) ? 20 : 50;
00085 ctrl.CoarsenTo = amin(vtxdist[npes]+1, vtx_factor*incon*amax(npes, inparts));
00086 ctrl.ipc_factor = iipc2redist;
00087 ctrl.redist_factor = 1.0;
00088 ctrl.redist_base = 1.0;
00089 ctrl.seed = (seed == 0 ? mype : seed*mype);
00090 ctrl.sync = GlobalSEMax(&ctrl, seed);
00091 ctrl.partType = ADAPTIVE_PARTITION;
00092 ctrl.ps_relation = ps_relation;
00093 ctrl.tpwgts = itpwgts;
00094
00095 graph = Moc_SetUpGraph(&ctrl, incon, vtxdist, xadj, vwgt, adjncy, adjwgt, &iwgtflag);
00096 graph->vsize = (vsize == NULL ? idxsmalloc(graph->nvtxs, 1, "vsize") : vsize);
00097
00098 graph->home = idxmalloc(graph->nvtxs, "home");
00099 if (ctrl.ps_relation == COUPLED)
00100 idxset(graph->nvtxs, mype, graph->home);
00101 else {
00102
00103 for (i=0; i<graph->nvtxs; i++)
00104 part[i] = (part[i] >= ctrl.nparts ? 0 : part[i]);
00105
00106 idxcopy(graph->nvtxs, part, graph->home);
00107 }
00108
00109 tewgt = idxsum(graph->nedges, graph->adjwgt);
00110 tvsize = idxsum(graph->nvtxs, graph->vsize);
00111 gtewgt = (floattype) GlobalSESum(&ctrl, tewgt) + 1.0/graph->gnvtxs;
00112 gtvsize = (floattype) GlobalSESum(&ctrl, tvsize) + 1.0/graph->gnvtxs;
00113 ctrl.edge_size_ratio = gtewgt/gtvsize;
00114 scopy(incon, iubvec, ctrl.ubvec);
00115
00116 PreAllocateMemory(&ctrl, graph, &wspace);
00117
00118
00119
00120
00121 IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl));
00122 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00123 IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr));
00124
00125 Adaptive_Partition(&ctrl, graph, &wspace);
00126 ParallelReMapGraph(&ctrl, graph, &wspace);
00127
00128 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00129 IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr));
00130
00131 idxcopy(graph->nvtxs, graph->where, part);
00132 if (edgecut != NULL)
00133 *edgecut = graph->mincut;
00134
00135
00136
00137
00138 IFSET(ctrl.dbglvl, DBG_TIME, PrintTimingInfo(&ctrl));
00139 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00140
00141 if (ctrl.dbglvl&DBG_INFO) {
00142 Mc_ComputeMoveStatistics(&ctrl, graph, &nmoved, &maxin, &maxout);
00143 rprintf(&ctrl, "Final %3d-way Cut: %6d \tBalance: ", inparts, graph->mincut);
00144 avg = 0.0;
00145 for (h=0; h<incon; h++) {
00146 maximb = 0.0;
00147 for (i=0; i<inparts; i++)
00148 maximb = amax(maximb, graph->gnpwgts[i*incon+h]/itpwgts[i*incon+h]);
00149 avg += maximb;
00150 rprintf(&ctrl, "%.3f ", maximb);
00151 }
00152 rprintf(&ctrl, "\nNMoved: %d %d %d %d\n", nmoved, maxin, maxout, maxin+maxout);
00153 }
00154
00155
00156
00157
00158 GKfree((void **)&graph->lnpwgts, (void **)&graph->gnpwgts, (void **)&graph->nvwgt, (void **)(&graph->home), LTERM);
00159 if (vsize == NULL)
00160 GKfree((void **)(&graph->vsize), LTERM);
00161 GKfree((void **)&itpwgts, LTERM);
00162 FreeInitialGraphAndRemap(graph, iwgtflag);
00163 FreeWSpace(&wspace);
00164 FreeCtrl(&ctrl);
00165
00166 if (inumflag == 1)
00167 ChangeNumbering(vtxdist, xadj, adjncy, part, npes, mype, 0);
00168
00169 return;
00170 }
00171
00172
00173
00174
00175
00176
00177
00178 void Adaptive_Partition(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00179 {
00180 int i;
00181 int tewgt, tvsize;
00182 floattype gtewgt, gtvsize;
00183 floattype ubavg, lbavg, lbvec[MAXNCON];
00184
00185
00186
00187
00188 SetUp(ctrl, graph, wspace);
00189
00190 ubavg = savg(graph->ncon, ctrl->ubvec);
00191 tewgt = idxsum(graph->nedges, graph->adjwgt);
00192 tvsize = idxsum(graph->nvtxs, graph->vsize);
00193 gtewgt = (floattype) GlobalSESum(ctrl, tewgt) + 1.0/graph->gnvtxs;
00194 gtvsize = (floattype) GlobalSESum(ctrl, tvsize) + 1.0/graph->gnvtxs;
00195 ctrl->redist_factor = ctrl->redist_base * ((gtewgt/gtvsize)/ ctrl->edge_size_ratio);
00196
00197 IFSET(ctrl->dbglvl, DBG_PROGRESS, rprintf(ctrl, "[%6d %8d %5d %5d][%d]\n",
00198 graph->gnvtxs, GlobalSESum(ctrl, graph->nedges), GlobalSEMin(ctrl, graph->nvtxs), GlobalSEMax(ctrl, graph->nvtxs), ctrl->CoarsenTo));
00199
00200 if (graph->gnvtxs < 1.3*ctrl->CoarsenTo ||
00201 (graph->finer != NULL && graph->gnvtxs > graph->finer->gnvtxs*COARSEN_FRACTION)) {
00202
00203
00204
00205
00206 graph->where = idxsmalloc(graph->nvtxs+graph->nrecv, -1, "graph->where");
00207 idxcopy(graph->nvtxs, graph->home, graph->where);
00208
00209 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00210 lbavg = savg(graph->ncon, lbvec);
00211
00212 if (lbavg > ubavg + 0.035 && ctrl->partType != REFINE_PARTITION)
00213 Balance_Partition(ctrl, graph, wspace);
00214
00215 if (ctrl->dbglvl&DBG_PROGRESS) {
00216 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00217 rprintf(ctrl, "nvtxs: %10d, balance: ", graph->gnvtxs);
00218 for (i=0; i<graph->ncon; i++)
00219 rprintf(ctrl, "%.3f ", lbvec[i]);
00220 rprintf(ctrl, "\n");
00221 }
00222
00223
00224 if (graph->finer == NULL) {
00225 Moc_ComputePartitionParams(ctrl, graph, wspace);
00226 Moc_KWayBalance(ctrl, graph, wspace, graph->ncon);
00227 Moc_KWayAdaptiveRefine(ctrl, graph, wspace, NGR_PASSES);
00228 }
00229 }
00230 else {
00231
00232
00233
00234 switch (ctrl->ps_relation) {
00235 case COUPLED:
00236 Mc_LocalMatch_HEM(ctrl, graph, wspace);
00237 break;
00238 case DISCOUPLED:
00239 default:
00240 Moc_GlobalMatch_Balance(ctrl, graph, wspace);
00241 break;
00242 }
00243
00244 Adaptive_Partition(ctrl, graph->coarser, wspace);
00245
00246
00247
00248
00249 Moc_ProjectPartition(ctrl, graph, wspace);
00250 Moc_ComputePartitionParams(ctrl, graph, wspace);
00251
00252 if (graph->ncon > 1 && graph->level < 4) {
00253 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00254 lbavg = savg(graph->ncon, lbvec);
00255
00256 if (lbavg > ubavg + 0.025) {
00257 Moc_KWayBalance(ctrl, graph, wspace, graph->ncon);
00258 }
00259 }
00260
00261 Moc_KWayAdaptiveRefine(ctrl, graph, wspace, NGR_PASSES);
00262
00263 if (ctrl->dbglvl&DBG_PROGRESS) {
00264 Moc_ComputeParallelBalance(ctrl, graph, graph->where, lbvec);
00265 rprintf(ctrl, "nvtxs: %10d, cut: %8d, balance: ", graph->gnvtxs, graph->mincut);
00266 for (i=0; i<graph->ncon; i++)
00267 rprintf(ctrl, "%.3f ", lbvec[i]);
00268 rprintf(ctrl, "\n");
00269 }
00270 }
00271 }
00272