00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #define DEBUG_ORDER_
00016
00017 #include <parmetislib.h>
00018
00019
00020
00021
00022 void MultilevelOrder(CtrlType *ctrl, GraphType *graph, idxtype *order, idxtype *sizes, WorkSpaceType *wspace)
00023 {
00024 int i, nparts, nvtxs, npes;
00025 idxtype *perm, *lastnode, *morder, *porder;
00026 GraphType *mgraph;
00027
00028 npes = ctrl->npes;
00029 nvtxs = graph->nvtxs;
00030
00031 perm = idxmalloc(nvtxs, "MultilevelOrder: perm");
00032 lastnode = idxsmalloc(4*npes, -1, "MultilevelOrder: lastnode");
00033
00034 for (i=0; i<nvtxs; i++)
00035 perm[i] = i;
00036 lastnode[2] = graph->gnvtxs;
00037
00038 idxset(nvtxs, -1, order);
00039
00040 sizes[0] = 2*npes-1;
00041
00042 graph->where = idxsmalloc(nvtxs, 0, "MultilevelOrder: graph->where");
00043
00044 for (nparts=2; nparts<=ctrl->npes; nparts*=2) {
00045 ctrl->nparts = nparts;
00046
00047 Order_Partition(ctrl, graph, wspace);
00048
00049 LabelSeparators(ctrl, graph, lastnode, perm, order, sizes, wspace);
00050
00051 CompactGraph(ctrl, graph, perm, wspace);
00052
00053 if (ctrl->CoarsenTo < 100*nparts) {
00054 ctrl->CoarsenTo = 1.5*ctrl->CoarsenTo;
00055 }
00056 ctrl->CoarsenTo = amin(ctrl->CoarsenTo, graph->gnvtxs-1);
00057 }
00058
00059
00060
00061
00062
00063 IFSET(ctrl->dbglvl, DBG_TIME, MPI_Barrier(ctrl->comm));
00064 IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->MoveTmr));
00065
00066 SetUp(ctrl, graph, wspace);
00067 graph->ncon = 1;
00068 mgraph = Moc_MoveGraph(ctrl, graph, wspace);
00069
00070
00071 for (i=0; i<npes; i++)
00072 sizes[i] = mgraph->vtxdist[i+1]-mgraph->vtxdist[i];
00073
00074 porder = idxmalloc(graph->nvtxs, "MultilevelOrder: porder");
00075 morder = idxmalloc(mgraph->nvtxs, "MultilevelOrder: morder");
00076
00077 IFSET(ctrl->dbglvl, DBG_TIME, MPI_Barrier(ctrl->comm));
00078 IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->MoveTmr));
00079
00080
00081 LocalNDOrder(ctrl, mgraph, morder, lastnode[2*(ctrl->npes+ctrl->mype)]-mgraph->nvtxs, wspace);
00082
00083
00084 ProjectInfoBack(ctrl, graph, porder, morder, wspace);
00085
00086
00087 for (i=0; i<graph->nvtxs; i++) {
00088 ASSERT(ctrl, order[perm[i]] == -1);
00089 order[perm[i]] = porder[i];
00090 }
00091
00092 FreeGraph(mgraph);
00093 GKfree((void **)&perm, (void **)&lastnode, (void **)&porder, (void **)&morder, LTERM);
00094
00095
00096 }
00097
00098
00099
00100
00101
00102
00103
00104 void LabelSeparators(CtrlType *ctrl, GraphType *graph, idxtype *lastnode, idxtype *perm, idxtype *order, idxtype *sizes, WorkSpaceType *wspace)
00105 {
00106 int i, nvtxs, nparts, sid;
00107 idxtype *where, *lpwgts, *gpwgts, *sizescan;
00108
00109 nparts = ctrl->nparts;
00110
00111 nvtxs = graph->nvtxs;
00112 where = graph->where;
00113 lpwgts = graph->lpwgts;
00114 gpwgts = graph->gpwgts;
00115
00116
00117
00118 idxset(2*nparts, 0, lpwgts);
00119 for (i=0; i<nvtxs; i++)
00120 lpwgts[where[i]]++;
00121
00122 sizescan = idxmalloc(2*nparts, "LabelSeparators: sizescan");
00123
00124
00125 MPI_Scan((void *)lpwgts, (void *)sizescan, 2*nparts, IDX_DATATYPE, MPI_SUM, ctrl->comm);
00126 MPI_Allreduce((void *)lpwgts, (void *)gpwgts, 2*nparts, IDX_DATATYPE, MPI_SUM, ctrl->comm);
00127
00128 #ifdef DEBUG_ORDER
00129 PrintVector(ctrl, 2*nparts, 0, lpwgts, "Lpwgts");
00130 PrintVector(ctrl, 2*nparts, 0, sizescan, "SizeScan");
00131 PrintVector(ctrl, 2*nparts, 0, lastnode, "LastNode");
00132 #endif
00133
00134
00135 for (i=nparts-2; i>=0; i-=2)
00136 sizes[--sizes[0]] = gpwgts[nparts+i];
00137
00138 if (ctrl->dbglvl&DBG_INFO) {
00139 if (ctrl->mype == 0) {
00140 printf("SepSizes: ");
00141 for (i=0; i<nparts; i+=2)
00142 printf(" %d [%d %d]", gpwgts[nparts+i], gpwgts[i], gpwgts[i+1]);
00143 printf("\n");
00144 }
00145 MPI_Barrier(ctrl->comm);
00146 }
00147
00148 for (i=0; i<2*nparts; i++)
00149 sizescan[i] -= lpwgts[i];
00150
00151 for (i=0; i<nvtxs; i++) {
00152 if (where[i] >= nparts) {
00153 sid = where[i];
00154 sizescan[sid]++;
00155 ASSERT(ctrl, order[perm[i]] == -1);
00156 order[perm[i]] = lastnode[sid] - sizescan[sid];
00157
00158 }
00159 }
00160
00161
00162 idxcopy(2*nparts, lastnode, sizescan);
00163 for (i=0; i<nparts; i+=2) {
00164 lastnode[2*nparts+2*i] = sizescan[nparts+i]-gpwgts[nparts+i]-gpwgts[i+1];
00165 lastnode[2*nparts+2*(i+1)] = sizescan[nparts+i]-gpwgts[nparts+i];
00166 }
00167
00168 free(sizescan);
00169
00170 }
00171
00172
00173
00174
00175
00176
00177
00178 void CompactGraph(CtrlType *ctrl, GraphType *graph, idxtype *perm, WorkSpaceType *wspace)
00179 {
00180 int i, j, l, nvtxs, cnvtxs, cfirstvtx, nparts, npes;
00181 idxtype *xadj, *ladjncy, *adjwgt, *vtxdist, *where;
00182 idxtype *cmap, *cvtxdist, *newwhere;
00183
00184 nparts = ctrl->nparts;
00185 npes = ctrl->npes;
00186
00187 nvtxs = graph->nvtxs;
00188 xadj = graph->xadj;
00189 ladjncy = graph->adjncy;
00190 adjwgt = graph->adjwgt;
00191 where = graph->where;
00192
00193 if (graph->cmap == NULL)
00194 graph->cmap = idxmalloc(nvtxs+graph->nrecv, "CompactGraph: cmap");
00195 cmap = graph->cmap;
00196
00197 vtxdist = graph->vtxdist;
00198
00199
00200
00201
00202
00203 cvtxdist = wspace->pv1;
00204 cnvtxs = cvtxdist[npes] = idxsum(nparts, graph->lpwgts);
00205
00206 MPI_Allgather((void *)(cvtxdist+npes), 1, IDX_DATATYPE, (void *)cvtxdist, 1, IDX_DATATYPE, ctrl->comm);
00207 MAKECSR(i, npes, cvtxdist);
00208
00209 #ifdef DEBUG_ORDER
00210 PrintVector(ctrl, npes+1, 0, cvtxdist, "cvtxdist");
00211 #endif
00212
00213
00214
00215
00216
00217 cfirstvtx = cvtxdist[ctrl->mype];
00218
00219
00220 for (cnvtxs=0, i=0; i<nvtxs; i++) {
00221 if (where[i] < nparts) {
00222 perm[cnvtxs] = perm[i];
00223 cmap[i] = cfirstvtx + cnvtxs++;
00224 }
00225 }
00226
00227 CommInterfaceData(ctrl, graph, cmap, wspace->indices, cmap+nvtxs);
00228
00229
00230
00231
00232
00233 newwhere = idxmalloc(cnvtxs, "CompactGraph: newwhere");
00234 cnvtxs = l = 0;
00235 for (i=0; i<nvtxs; i++) {
00236 if (where[i] < nparts) {
00237 for (j=xadj[i]; j<xadj[i+1]; j++) {
00238 if (where[i] == where[ladjncy[j]]) {
00239 ladjncy[l] = cmap[ladjncy[j]];
00240 adjwgt[l++] = adjwgt[j];
00241 }
00242 #ifdef DEBUG_ORDER
00243 else if (where[ladjncy[j]] < nparts)
00244 printf("It seems that the separation has failed: %d %d\n", where[i], where[ladjncy[j]]);
00245 #endif
00246 }
00247
00248 xadj[cnvtxs] = l;
00249 graph->vwgt[cnvtxs] = graph->vwgt[i];
00250 newwhere[cnvtxs] = where[i];
00251 cnvtxs++;
00252 }
00253 }
00254 for (i=cnvtxs; i>0; i--)
00255 xadj[i] = xadj[i-1];
00256 xadj[0] = 0;
00257
00258 GKfree((void **)&graph->match, (void **)&graph->cmap, (void **)&graph->lperm, (void **)&graph->where, (void **)&graph->label, (void **)&graph->rinfo,
00259 (void **)&graph->nrinfo, (void **)&graph->lpwgts, (void **)&graph->gpwgts, (void **)&graph->sepind, (void **)&graph->peind,
00260 (void **)&graph->sendptr, (void **)&graph->sendind, (void **)&graph->recvptr, (void **)&graph->recvind,
00261 (void **)&graph->imap, (void **)&graph->rlens, (void **)&graph->slens, (void **)&graph->rcand, (void **)&graph->pexadj,
00262 (void **)&graph->peadjncy, (void **)&graph->peadjloc, LTERM);
00263
00264 graph->nvtxs = cnvtxs;
00265 graph->nedges = l;
00266 graph->gnvtxs = cvtxdist[npes];
00267 idxcopy(npes+1, cvtxdist, graph->vtxdist);
00268 graph->where = newwhere;
00269
00270 }
00271
00272
00273
00274
00275
00276
00277 void LocalNDOrder(CtrlType *ctrl, GraphType *graph, idxtype *order, int firstnode, WorkSpaceType *wspace)
00278 {
00279 int i, j, nvtxs, firstvtx, lastvtx;
00280 idxtype *xadj, *adjncy;
00281 idxtype *perm, *iperm;
00282 int numflag=0, options[10];
00283
00284 nvtxs = graph->nvtxs;
00285 xadj = graph->xadj;
00286 adjncy = graph->adjncy;
00287
00288 firstvtx = graph->vtxdist[ctrl->mype];
00289 lastvtx = graph->vtxdist[ctrl->mype+1];
00290
00291
00292 for (i=0; i<nvtxs; i++) {
00293 for (j=xadj[i]; j<xadj[i+1]; j++) {
00294 ASSERT(ctrl, adjncy[j]>=firstvtx && adjncy[j]<lastvtx);
00295 adjncy[j] -= firstvtx;
00296 }
00297 }
00298
00299 ASSERT(ctrl, 2*(nvtxs+5) < wspace->maxcore);
00300
00301 perm = wspace->core;
00302 iperm = perm + nvtxs + 5;
00303
00304 options[0] = 0;
00305 METIS_NodeND(&nvtxs, xadj, adjncy, &numflag, options, perm, iperm);
00306
00307 for (i=0; i<nvtxs; i++) {
00308 ASSERT(ctrl, iperm[i]>=0 && iperm[i]<nvtxs);
00309 order[i] = firstnode+iperm[i];
00310 }
00311
00312 }
00313
00314
00315
00316
00317 void Order_Partition(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00318 {
00319
00320 SetUp(ctrl, graph, wspace);
00321 graph->ncon = 1;
00322
00323 IFSET(ctrl->dbglvl, DBG_PROGRESS, rprintf(ctrl, "[%6d %8d %5d %5d][%d][%d]\n",
00324 graph->gnvtxs, GlobalSESum(ctrl, graph->nedges), GlobalSEMin(ctrl, graph->nvtxs),
00325 GlobalSEMax(ctrl, graph->nvtxs), ctrl->CoarsenTo,
00326 GlobalSEMax(ctrl, graph->vwgt[idxamax(graph->nvtxs, graph->vwgt)])));
00327
00328 if (graph->gnvtxs < 1.3*ctrl->CoarsenTo || (graph->finer != NULL && graph->gnvtxs > graph->finer->gnvtxs*COARSEN_FRACTION)) {
00329
00330 InitMultisection(ctrl, graph, wspace);
00331
00332 if (graph->finer == NULL) {
00333 ComputeNodePartitionParams(ctrl, graph, wspace);
00334 KWayNodeRefine(ctrl, graph, wspace, 2*NGR_PASSES, ORDER_UNBALANCE_FRACTION);
00335 }
00336 }
00337 else {
00338 Mc_LocalMatch_HEM(ctrl, graph, wspace);
00339
00340 Order_Partition(ctrl, graph->coarser, wspace);
00341
00342 Moc_ProjectPartition(ctrl, graph, wspace);
00343 ComputeNodePartitionParams(ctrl, graph, wspace);
00344 KWayNodeRefine(ctrl, graph, wspace, 2*NGR_PASSES, ORDER_UNBALANCE_FRACTION);
00345 }
00346 }
00347
00348