/************************************************************************* * This function is the entry point for PWMETIS that accepts exact weights * for the target partitions **************************************************************************/ void METIS_WPartGraphRecursive(int *nvtxs, idxtype *xadj, idxtype *adjncy, idxtype *vwgt, idxtype *adjwgt, int *wgtflag, int *numflag, int *nparts, floattype *tpwgts, int *options, int *edgecut, idxtype *part) { int i, j; GraphType graph; CtrlType ctrl; floattype *mytpwgts; if (*numflag == 1) Change2CNumbering(*nvtxs, xadj, adjncy); SetUpGraph(&graph, OP_PMETIS, *nvtxs, 1, xadj, adjncy, vwgt, adjwgt, *wgtflag); if (options[0] == 0) { /* Use the default parameters */ ctrl.CType = PMETIS_CTYPE; ctrl.IType = PMETIS_ITYPE; ctrl.RType = PMETIS_RTYPE; ctrl.dbglvl = PMETIS_DBGLVL; } else { ctrl.CType = options[OPTION_CTYPE]; ctrl.IType = options[OPTION_ITYPE]; ctrl.RType = options[OPTION_RTYPE]; ctrl.dbglvl = options[OPTION_DBGLVL]; } ctrl.optype = OP_PMETIS; ctrl.CoarsenTo = 20; ctrl.maxvwgt = 1.5*(idxsum(*nvtxs, graph.vwgt)/ctrl.CoarsenTo); mytpwgts = fmalloc(*nparts, "PWMETIS: mytpwgts"); for (i=0; i<*nparts; i++) mytpwgts[i] = tpwgts[i]; InitRandom(-1); AllocateWorkSpace(&ctrl, &graph, *nparts); IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl)); IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr)); *edgecut = MlevelRecursiveBisection(&ctrl, &graph, *nparts, part, mytpwgts, 1.000, 0); IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr)); IFSET(ctrl.dbglvl, DBG_TIME, PrintTimers(&ctrl)); FreeWorkSpace(&ctrl, &graph); free(mytpwgts); if (*numflag == 1) Change2FNumbering(*nvtxs, xadj, adjncy, part); }
/************************************************************************* * This function is the entry point for KWMETIS **************************************************************************/ void METIS_mCPartGraphKway(int *nvtxs, int *ncon, idxtype *xadj, idxtype *adjncy, idxtype *vwgt, idxtype *adjwgt, int *wgtflag, int *numflag, int *nparts, floattype *rubvec, int *options, int *edgecut, idxtype *part) { int i, j; GraphType graph; CtrlType ctrl; if (*numflag == 1) Change2CNumbering(*nvtxs, xadj, adjncy); SetUpGraph(&graph, OP_KMETIS, *nvtxs, *ncon, xadj, adjncy, vwgt, adjwgt, *wgtflag); if (options[0] == 0) { /* Use the default parameters */ ctrl.CType = McKMETIS_CTYPE; ctrl.IType = McKMETIS_ITYPE; ctrl.RType = McKMETIS_RTYPE; ctrl.dbglvl = McKMETIS_DBGLVL; } else { ctrl.CType = options[OPTION_CTYPE]; ctrl.IType = options[OPTION_ITYPE]; ctrl.RType = options[OPTION_RTYPE]; ctrl.dbglvl = options[OPTION_DBGLVL]; } ctrl.optype = OP_KMETIS; ctrl.CoarsenTo = amax((*nvtxs)/(20*log2Int(*nparts)), 30*(*nparts)); ctrl.nmaxvwgt = 1.5/(1.0*ctrl.CoarsenTo); InitRandom(-1); AllocateWorkSpace(&ctrl, &graph, *nparts); IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl)); IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr)); ASSERT(CheckGraph(&graph)); *edgecut = MCMlevelKWayPartitioning(&ctrl, &graph, *nparts, part, rubvec); IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr)); IFSET(ctrl.dbglvl, DBG_TIME, PrintTimers(&ctrl)); FreeWorkSpace(&ctrl, &graph); if (*numflag == 1) Change2FNumbering(*nvtxs, xadj, adjncy, part); }
/************************************************************************* * This function is the entry point for KWMETIS **************************************************************************/ void METIS_WPartGraphVKway(int *nvtxs, idxtype *xadj, idxtype *adjncy, idxtype *vwgt, idxtype *vsize, int *wgtflag, int *numflag, int *nparts, float *tpwgts, int *options, int *volume, idxtype *part) { int i, j; GraphType graph; CtrlType ctrl; if (*numflag == 1) Change2CNumbering(*nvtxs, xadj, adjncy); VolSetUpGraph(&graph, OP_KVMETIS, *nvtxs, 1, xadj, adjncy, vwgt, vsize, *wgtflag); if (options[0] == 0) { /* Use the default parameters */ ctrl.CType = KVMETIS_CTYPE; ctrl.IType = KVMETIS_ITYPE; ctrl.RType = KVMETIS_RTYPE; ctrl.dbglvl = KVMETIS_DBGLVL; } else { ctrl.CType = options[OPTION_CTYPE]; ctrl.IType = options[OPTION_ITYPE]; ctrl.RType = options[OPTION_RTYPE]; ctrl.dbglvl = options[OPTION_DBGLVL]; } ctrl.optype = OP_KVMETIS; ctrl.CoarsenTo = amax((*nvtxs)/(40*log2Int(*nparts)), 20*(*nparts)); ctrl.maxvwgt = 1.5*((graph.vwgt ? idxsum(*nvtxs, graph.vwgt) : (*nvtxs))/ctrl.CoarsenTo); InitRandom(-1); AllocateWorkSpace(&ctrl, &graph, *nparts); IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl)); IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr)); *volume = MlevelVolKWayPartitioning(&ctrl, &graph, *nparts, part, tpwgts, 1.03); IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr)); IFSET(ctrl.dbglvl, DBG_TIME, PrintTimers(&ctrl)); FreeWorkSpace(&ctrl, &graph); if (*numflag == 1) Change2FNumbering(*nvtxs, xadj, adjncy, part); }
/************************************************************************* * This function is the entry point for KWMETIS with seed specification * in options[7] **************************************************************************/ void METIS_WPartGraphKway2(idxtype *nvtxs, idxtype *xadj, idxtype *adjncy, idxtype *vwgt, idxtype *adjwgt, idxtype *wgtflag, idxtype *numflag, idxtype *nparts, float *tpwgts, idxtype *options, idxtype *edgecut, idxtype *part) { idxtype i, j; GraphType graph; CtrlType ctrl; if (*numflag == 1) Change2CNumbering(*nvtxs, xadj, adjncy); SetUpGraph(&graph, OP_KMETIS, *nvtxs, 1, xadj, adjncy, vwgt, adjwgt, *wgtflag); if (options[0] == 0) { /* Use the default parameters */ ctrl.CType = KMETIS_CTYPE; ctrl.IType = KMETIS_ITYPE; ctrl.RType = KMETIS_RTYPE; ctrl.dbglvl = KMETIS_DBGLVL; } else { ctrl.CType = options[OPTION_CTYPE]; ctrl.IType = options[OPTION_ITYPE]; ctrl.RType = options[OPTION_RTYPE]; ctrl.dbglvl = options[OPTION_DBGLVL]; } ctrl.optype = OP_KMETIS; ctrl.CoarsenTo = 20*(*nparts); ctrl.maxvwgt = 1.5*((graph.vwgt ? idxsum(*nvtxs, graph.vwgt, 1) : (*nvtxs))/ctrl.CoarsenTo); InitRandom(options[7]); AllocateWorkSpace(&ctrl, &graph, *nparts); IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl)); IFSET(ctrl.dbglvl, DBG_TIME, gk_startcputimer(ctrl.TotalTmr)); *edgecut = MlevelKWayPartitioning(&ctrl, &graph, *nparts, part, tpwgts, 1.03); IFSET(ctrl.dbglvl, DBG_TIME, gk_stopcputimer(ctrl.TotalTmr)); IFSET(ctrl.dbglvl, DBG_TIME, PrintTimers(&ctrl)); FreeWorkSpace(&ctrl, &graph); if (*numflag == 1) Change2FNumbering(*nvtxs, xadj, adjncy, part); }
int METIS_PartGraphKway(idx_t *nvtxs, idx_t *ncon, idx_t *xadj, idx_t *adjncy, idx_t *vwgt, idx_t *vsize, idx_t *adjwgt, idx_t *nparts, real_t *tpwgts, real_t *ubvec, idx_t *options, idx_t *objval, idx_t *part) { int sigrval=0, renumber=0; graph_t *graph; ctrl_t *ctrl; /* set up malloc cleaning code and signal catchers */ if (!gk_malloc_init()) return METIS_ERROR_MEMORY; gk_sigtrap(); if ((sigrval = gk_sigcatch()) != 0) goto SIGTHROW; /* set up the run parameters */ ctrl = SetupCtrl(METIS_OP_KMETIS, options, *ncon, *nparts, tpwgts, ubvec); if (!ctrl) { gk_siguntrap(); return METIS_ERROR_INPUT; } /* if required, change the numbering to 0 */ if (ctrl->numflag == 1) { Change2CNumbering(*nvtxs, xadj, adjncy); renumber = 1; } /* set up the graph */ graph = SetupGraph(ctrl, *nvtxs, *ncon, xadj, adjncy, vwgt, vsize, adjwgt); /* set up multipliers for making balance computations easier */ SetupKWayBalMultipliers(ctrl, graph); /* set various run parameters that depend on the graph */ if (ctrl->iptype == METIS_IPTYPE_METISRB) { ctrl->CoarsenTo = gk_max((*nvtxs)/(40*gk_log2(*nparts)), 30*(*nparts)); ctrl->CoarsenTo = 10*(*nparts); ctrl->nIparts = (ctrl->CoarsenTo == 30*(*nparts) ? 4 : 5); } else { ctrl->CoarsenTo = 10*(*nparts); ctrl->nIparts = 10; } /* take care contiguity requests for disconnected graphs */ if (ctrl->contig && !IsConnected(graph, 0)) gk_errexit(SIGERR, "METIS Error: A contiguous partition is requested for a non-contiguous input graph.\n"); /* allocate workspace memory */ AllocateWorkSpace(ctrl, graph); /* start the partitioning */ IFSET(ctrl->dbglvl, METIS_DBG_TIME, InitTimers(ctrl)); IFSET(ctrl->dbglvl, METIS_DBG_TIME, gk_startwctimer(ctrl->TotalTmr)); *objval = MlevelKWayPartitioning(ctrl, graph, part); IFSET(ctrl->dbglvl, METIS_DBG_TIME, gk_stopwctimer(ctrl->TotalTmr)); IFSET(ctrl->dbglvl, METIS_DBG_TIME, PrintTimers(ctrl)); /* clean up */ FreeCtrl(&ctrl); SIGTHROW: /* if required, change the numbering back to 1 */ if (renumber) Change2FNumbering(*nvtxs, xadj, adjncy, part); gk_siguntrap(); gk_malloc_cleanup(0); return metis_rcode(sigrval); }
/************************************************************************* * This function is the entry point for PWMETIS that accepts exact weights * for the target partitions **************************************************************************/ void METIS_mCPartGraphRecursive2(int *nvtxs, int *ncon, idxtype *xadj, idxtype *adjncy, idxtype *vwgt, idxtype *adjwgt, int *wgtflag, int *numflag, int *nparts, float *tpwgts, int *options, int *edgecut, idxtype *part) { int i, j; GraphType graph; CtrlType ctrl; float *mytpwgts; float avgwgt; if (*numflag == 1) Change2CNumbering(*nvtxs, xadj, adjncy); SetUpGraph(&graph, OP_PMETIS, *nvtxs, *ncon, xadj, adjncy, vwgt, adjwgt, *wgtflag); graph.npwgts = NULL; mytpwgts = fmalloc(*nparts, "mytpwgts"); scopy(*nparts, tpwgts, mytpwgts); if (options[0] == 0) { /* Use the default parameters */ ctrl.CType = McPMETIS_CTYPE; ctrl.IType = McPMETIS_ITYPE; ctrl.RType = McPMETIS_RTYPE; ctrl.dbglvl = McPMETIS_DBGLVL; } else { ctrl.CType = options[OPTION_CTYPE]; ctrl.IType = options[OPTION_ITYPE]; ctrl.RType = options[OPTION_RTYPE]; ctrl.dbglvl = options[OPTION_DBGLVL]; } ctrl.optype = OP_PMETIS; ctrl.CoarsenTo = 100; ctrl.nmaxvwgt = 1.5/(1.0*ctrl.CoarsenTo); InitRandom(options[7]); AllocateWorkSpace(&ctrl, &graph, *nparts); IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl)); IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr)); ASSERT(CheckGraph(&graph)); *edgecut = MCMlevelRecursiveBisection2(&ctrl, &graph, *nparts, mytpwgts, part, 1.000, 0); /* { idxtype wgt[2048], minwgt, maxwgt, sumwgt; printf("nvtxs: %d, nparts: %d, ncon: %d\n", graph.nvtxs, *nparts, *ncon); for (i=0; i<(*nparts)*(*ncon); i++) wgt[i] = 0; for (i=0; i<graph.nvtxs; i++) for (j=0; j<*ncon; j++) wgt[part[i]*(*ncon)+j] += vwgt[i*(*ncon)+j]; for (j=0; j<*ncon; j++) { minwgt = maxwgt = sumwgt = 0; for (i=0; i<(*nparts); i++) { minwgt = (wgt[i*(*ncon)+j] < wgt[minwgt*(*ncon)+j]) ? i : minwgt; maxwgt = (wgt[i*(*ncon)+j] > wgt[maxwgt*(*ncon)+j]) ? i : maxwgt; sumwgt += wgt[i*(*ncon)+j]; } avgwgt = (float)sumwgt / (float)*nparts; printf("min: %5d, max: %5d, avg: %5.2f, balance: %6.3f\n", wgt[minwgt*(*ncon)+j], wgt[maxwgt*(*ncon)+j], avgwgt, (float)wgt[maxwgt*(*ncon)+j] / avgwgt); } printf("\n"); } */ IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr)); IFSET(ctrl.dbglvl, DBG_TIME, PrintTimers(&ctrl)); FreeWorkSpace(&ctrl, &graph); GKfree((void**)(void *)&mytpwgts, LTERM); if (*numflag == 1) Change2FNumbering(*nvtxs, xadj, adjncy, part); }
/************************************************************************* * This function is the entry point for KMETIS **************************************************************************/ void *METIS_PartGraphForContact(idxtype *nvtxs, idxtype *xadj, idxtype *adjncy, double *xyzcoords, idxtype *sflag, idxtype *numflag, idxtype *nparts, idxtype *options, idxtype *edgecut, idxtype *part) { idxtype i, j, ii, dim, ncon, wgtflag, mcnumflag, nnodes, nlnodes, nclean, naclean, ndirty, maxdepth, rwgtflag, rnumflag; idxtype *mcvwgt, *dtpart, *marker, *leafpart; idxtype *adjwgt; float rubvec[2], lbvec[2]; GraphType graph, *cgraph; ContactInfoType *cinfo; DKeyValueType *xyzcand[3]; if (*numflag == 1) Change2CNumbering(*nvtxs, xadj, adjncy); /*--------------------------------------------------------------------- * Allocate memory for the contact info type *---------------------------------------------------------------------*/ cinfo = (ContactInfoType *)gk_malloc(sizeof(ContactInfoType), "METIS_PartGraphForContact: cinfo"); cinfo->leafptr = idxsmalloc(*nvtxs+1, 0, "METIS_PartGraphForContact: leafptr"); cinfo->leafind = idxsmalloc(*nvtxs, 0, "METIS_PartGraphForContact: leafind"); cinfo->leafwgt = idxsmalloc(*nvtxs, 0, "METIS_PartGraphForContact: leafwgt"); cinfo->part = idxsmalloc(*nvtxs, 0, "METIS_PartGraphForContact: part"); leafpart = cinfo->leafpart = idxmalloc(*nvtxs, "METIS_PartGraphForContact: leafpart"); cinfo->dtree = (DTreeNodeType *)gk_malloc(sizeof(DTreeNodeType)*(*nvtxs), "METIS_PartGraphForContact: cinfo->dtree"); cinfo->nvtxs = *nvtxs; /*--------------------------------------------------------------------- * Compute the initial k-way partitioning *---------------------------------------------------------------------*/ mcvwgt = idxsmalloc(2*(*nvtxs), 0, "METIS_PartGraphForContact: mcvwgt"); for (i=0; i<*nvtxs; i++) { mcvwgt[2*i+0] = 1; mcvwgt[2*i+1] = (sflag[i] == 0 ? 0 : 1); } adjwgt = idxmalloc(xadj[*nvtxs], "METIS_PartGraphForContact: adjwgt"); for (i=0; i<*nvtxs; i++) { for (j=xadj[i]; j<xadj[i+1]; j++) adjwgt[j] = (sflag[i] && sflag[adjncy[j]] ? 5 : 1); } rubvec[0] = 1.03; rubvec[1] = 1.05; ncon = 2; mcnumflag = 0; wgtflag = 1; METIS_mCPartGraphKway(nvtxs, &ncon, xadj, adjncy, mcvwgt, adjwgt, &wgtflag, &mcnumflag, nparts, rubvec, options, edgecut, part); /* The following is just for stat reporting purposes */ SetUpGraph(&graph, OP_KMETIS, *nvtxs, 2, xadj, adjncy, mcvwgt, NULL, 0); graph.vwgt = mcvwgt; ComputePartitionBalance(&graph, *nparts, part, lbvec); mprintf(" %D-way Edge-Cut: %7D, Balance: %5.2f %5.2f\n", *nparts, ComputeCut(&graph, part), lbvec[0], lbvec[1]); /*--------------------------------------------------------------------- * Induce the decission tree *---------------------------------------------------------------------*/ dtpart = idxmalloc(*nvtxs, "METIS_PartGraphForContact: dtpart"); marker = idxsmalloc(*nvtxs, 0, "METIS_PartGraphForContact: marker"); for (dim=0; dim<3; dim++) { xyzcand[dim] = (DKeyValueType *)gk_malloc(sizeof(DKeyValueType)*(*nvtxs), "METIS_PartGraphForContact: xyzcand[dim]"); for (i=0; i<*nvtxs; i++) { xyzcand[dim][i].key = xyzcoords[3*i+dim]; xyzcand[dim][i].val = i; } idkeysort(*nvtxs, xyzcand[dim]); } nnodes = nlnodes = nclean = naclean = ndirty = maxdepth = 0; InduceDecissionTree(*nvtxs, xyzcand, sflag, *nparts, part, *nvtxs/(20*(*nparts)), *nvtxs/(20*(*nparts)*(*nparts)), 0.90, &nnodes, &nlnodes, cinfo->dtree, leafpart, dtpart, &nclean, &naclean, &ndirty, &maxdepth, marker); mprintf("NNodes: %5D, NLNodes: %5D, NClean: %5D, NAClean: %5D, NDirty: %5D, MaxDepth: %3D\n", nnodes, nlnodes, nclean, naclean, ndirty, maxdepth); /*--------------------------------------------------------------------- * Create the tree-induced coarse graph and refine it *---------------------------------------------------------------------*/ cgraph = CreatePartitionGraphForContact(*nvtxs, xadj, adjncy, mcvwgt, adjwgt, nlnodes, leafpart); for (i=0; i<*nvtxs; i++) part[leafpart[i]] = dtpart[i]; ComputePartitionBalance(cgraph, *nparts, part, lbvec); mprintf(" %D-way Edge-Cut: %7D, Balance: %5.2f %5.2f\n", *nparts, ComputeCut(cgraph, part), lbvec[0], lbvec[1]); rwgtflag = 3; rnumflag = 0; METIS_mCRefineGraphKway(&(cgraph->nvtxs), &ncon, cgraph->xadj, cgraph->adjncy, cgraph->vwgt, cgraph->adjwgt, &rwgtflag, &rnumflag, nparts, rubvec, options, edgecut, part); ComputePartitionBalance(cgraph, *nparts, part, lbvec); mprintf(" %D-way Edge-Cut: %7D, Balance: %5.2f %5.2f\n", *nparts, ComputeCut(cgraph, part), lbvec[0], lbvec[1]); /*--------------------------------------------------------------------- * Use that to compute the partition of the original graph *---------------------------------------------------------------------*/ idxcopy(cgraph->nvtxs, part, dtpart); for (i=0; i<*nvtxs; i++) part[i] = dtpart[leafpart[i]]; ComputePartitionBalance(&graph, *nparts, part, lbvec); idxset(*nvtxs, 1, graph.vwgt); mprintf(" %D-way Edge-Cut: %7D, Volume: %7D, Balance: %5.2f %5.2f\n", *nparts, ComputeCut(&graph, part), ComputeVolume(&graph, part), lbvec[0], lbvec[1]); /*--------------------------------------------------------------------- * Induce the final decission tree *---------------------------------------------------------------------*/ nnodes = nlnodes = nclean = naclean = ndirty = maxdepth = 0; InduceDecissionTree(*nvtxs, xyzcand, sflag, *nparts, part, *nvtxs/((40)*(*nparts)), 1, 1.00, &nnodes, &nlnodes, cinfo->dtree, leafpart, dtpart, &nclean, &naclean, &ndirty, &maxdepth, marker); mprintf("NNodes: %5D, NLNodes: %5D, NClean: %5D, NAClean: %5D, NDirty: %5D, MaxDepth: %3D\n", nnodes, nlnodes, nclean, naclean, ndirty, maxdepth); /*--------------------------------------------------------------------- * Populate the remaining fields of the cinfo data structure *---------------------------------------------------------------------*/ cinfo->nnodes = nnodes; cinfo->nleafs = nlnodes; idxcopy(*nvtxs, part, cinfo->part); BuildDTLeafContents(cinfo, sflag); CheckDTree(*nvtxs, xyzcoords, part, cinfo); gk_free((void **)&mcvwgt, &dtpart, &xyzcand[0], &xyzcand[1], &xyzcand[2], &marker, &adjwgt, LTERM); if (*numflag == 1) Change2FNumbering(*nvtxs, xadj, adjncy, part); return (void *)cinfo; }