00001
00011 #include "metislib.h"
00012
00013
00014
00016
00017 ctrl_t *SetupCtrl(moptype_et optype, idx_t *options, idx_t ncon, idx_t nparts,
00018 real_t *tpwgts, real_t *ubvec)
00019 {
00020 idx_t i, j;
00021 ctrl_t *ctrl;
00022
00023 ctrl = (ctrl_t *)gk_malloc(sizeof(ctrl_t), "SetupCtrl: ctrl");
00024
00025 memset((void *)ctrl, 0, sizeof(ctrl_t));
00026
00027 switch (optype) {
00028 case METIS_OP_PMETIS:
00029 ctrl->objtype = GETOPTION(options, METIS_OPTION_OBJTYPE, METIS_OBJTYPE_CUT);
00030 ctrl->rtype = METIS_RTYPE_FM;
00031 ctrl->ncuts = GETOPTION(options, METIS_OPTION_NCUTS, 1);
00032 ctrl->niter = GETOPTION(options, METIS_OPTION_NITER, 10);
00033
00034 if (ncon == 1) {
00035 ctrl->iptype = GETOPTION(options, METIS_OPTION_IPTYPE, METIS_IPTYPE_GROW);
00036 ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, PMETIS_DEFAULT_UFACTOR);
00037 ctrl->CoarsenTo = 20;
00038 }
00039 else {
00040 ctrl->iptype = GETOPTION(options, METIS_OPTION_IPTYPE, METIS_IPTYPE_RANDOM);
00041 ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, MCPMETIS_DEFAULT_UFACTOR);
00042 ctrl->CoarsenTo = 100;
00043 }
00044
00045 break;
00046
00047
00048 case METIS_OP_KMETIS:
00049 ctrl->objtype = GETOPTION(options, METIS_OPTION_OBJTYPE, METIS_OBJTYPE_CUT);
00050 ctrl->iptype = METIS_IPTYPE_METISRB;
00051 ctrl->rtype = METIS_RTYPE_GREEDY;
00052 ctrl->ncuts = GETOPTION(options, METIS_OPTION_NCUTS, 1);
00053 ctrl->niter = GETOPTION(options, METIS_OPTION_NITER, 10);
00054 ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, KMETIS_DEFAULT_UFACTOR);
00055 ctrl->minconn = GETOPTION(options, METIS_OPTION_MINCONN, 0);
00056 ctrl->contig = GETOPTION(options, METIS_OPTION_CONTIG, 0);
00057 break;
00058
00059
00060 case METIS_OP_OMETIS:
00061 ctrl->objtype = GETOPTION(options, METIS_OPTION_OBJTYPE, METIS_OBJTYPE_NODE);
00062 ctrl->rtype = GETOPTION(options, METIS_OPTION_RTYPE, METIS_RTYPE_SEP1SIDED);
00063 ctrl->iptype = GETOPTION(options, METIS_OPTION_IPTYPE, METIS_IPTYPE_EDGE);
00064 ctrl->nseps = GETOPTION(options, METIS_OPTION_NSEPS, 1);
00065 ctrl->niter = GETOPTION(options, METIS_OPTION_NITER, 10);
00066 ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, OMETIS_DEFAULT_UFACTOR);
00067 ctrl->compress = GETOPTION(options, METIS_OPTION_COMPRESS, 1);
00068 ctrl->ccorder = GETOPTION(options, METIS_OPTION_CCORDER, 0);
00069 ctrl->pfactor = 0.1*GETOPTION(options, METIS_OPTION_PFACTOR, 0);
00070
00071 ctrl->CoarsenTo = 100;
00072 break;
00073
00074 default:
00075 gk_errexit(SIGERR, "Unknown optype of %d\n", optype);
00076 }
00077
00078
00079 ctrl->ctype = GETOPTION(options, METIS_OPTION_CTYPE, METIS_CTYPE_SHEM);
00080 ctrl->no2hop = GETOPTION(options, METIS_OPTION_NO2HOP, 0);
00081 ctrl->seed = GETOPTION(options, METIS_OPTION_SEED, -1);
00082 ctrl->dbglvl = GETOPTION(options, METIS_OPTION_DBGLVL, 0);
00083 ctrl->numflag = GETOPTION(options, METIS_OPTION_NUMBERING, 0);
00084
00085
00086 ctrl->optype = optype;
00087 ctrl->ncon = ncon;
00088 ctrl->nparts = nparts;
00089 ctrl->maxvwgt = ismalloc(ncon, 0, "SetupCtrl: maxvwgt");
00090
00091
00092 if (ctrl->optype != METIS_OP_OMETIS) {
00093 ctrl->tpwgts = rmalloc(nparts*ncon, "SetupCtrl: ctrl->tpwgts");
00094 if (tpwgts) {
00095 rcopy(nparts*ncon, tpwgts, ctrl->tpwgts);
00096 }
00097 else {
00098 for (i=0; i<nparts; i++) {
00099 for (j=0; j<ncon; j++)
00100 ctrl->tpwgts[i*ncon+j] = 1.0/nparts;
00101 }
00102 }
00103 }
00104 else {
00105
00106
00107 ctrl->tpwgts = rsmalloc(2, .5, "SetupCtrl: ctrl->tpwgts");
00108 }
00109
00110
00111
00112 ctrl->ubfactors = rsmalloc(ctrl->ncon, I2RUBFACTOR(ctrl->ufactor), "SetupCtrl: ubfactors");
00113 if (ubvec)
00114 rcopy(ctrl->ncon, ubvec, ctrl->ubfactors);
00115 for (i=0; i<ctrl->ncon; i++)
00116 ctrl->ubfactors[i] += 0.0000499;
00117
00118
00119
00120
00121 ctrl->pijbm = rmalloc(nparts*ncon, "SetupCtrl: ctrl->pijbm");
00122
00123 InitRandom(ctrl->seed);
00124
00125 IFSET(ctrl->dbglvl, METIS_DBG_INFO, PrintCtrl(ctrl));
00126
00127 if (!CheckParams(ctrl)) {
00128 FreeCtrl(&ctrl);
00129 return NULL;
00130 }
00131 else {
00132 return ctrl;
00133 }
00134 }
00135
00136
00137
00139
00140 void SetupKWayBalMultipliers(ctrl_t *ctrl, graph_t *graph)
00141 {
00142 idx_t i, j;
00143
00144 for (i=0; i<ctrl->nparts; i++) {
00145 for (j=0; j<graph->ncon; j++)
00146 ctrl->pijbm[i*graph->ncon+j] = graph->invtvwgt[j]/ctrl->tpwgts[i*graph->ncon+j];
00147 }
00148 }
00149
00150
00151
00153
00154 void Setup2WayBalMultipliers(ctrl_t *ctrl, graph_t *graph, real_t *tpwgts)
00155 {
00156 idx_t i, j;
00157
00158 for (i=0; i<2; i++) {
00159 for (j=0; j<graph->ncon; j++)
00160 ctrl->pijbm[i*graph->ncon+j] = graph->invtvwgt[j]/tpwgts[i*graph->ncon+j];
00161 }
00162 }
00163
00164
00165
00167
00168 void PrintCtrl(ctrl_t *ctrl)
00169 {
00170 idx_t i, j, modnum;
00171
00172 printf(" Runtime parameters:\n");
00173
00174 printf(" Objective type: ");
00175 switch (ctrl->objtype) {
00176 case METIS_OBJTYPE_CUT:
00177 printf("METIS_OBJTYPE_CUT\n");
00178 break;
00179 case METIS_OBJTYPE_VOL:
00180 printf("METIS_OBJTYPE_VOL\n");
00181 break;
00182 case METIS_OBJTYPE_NODE:
00183 printf("METIS_OBJTYPE_NODE\n");
00184 break;
00185 default:
00186 printf("Unknown!\n");
00187 }
00188
00189 printf(" Coarsening type: ");
00190 switch (ctrl->ctype) {
00191 case METIS_CTYPE_RM:
00192 printf("METIS_CTYPE_RM\n");
00193 break;
00194 case METIS_CTYPE_SHEM:
00195 printf("METIS_CTYPE_SHEM\n");
00196 break;
00197 default:
00198 printf("Unknown!\n");
00199 }
00200
00201 printf(" Initial partitioning type: ");
00202 switch (ctrl->iptype) {
00203 case METIS_IPTYPE_GROW:
00204 printf("METIS_IPTYPE_GROW\n");
00205 break;
00206 case METIS_IPTYPE_RANDOM:
00207 printf("METIS_IPTYPE_RANDOM\n");
00208 break;
00209 case METIS_IPTYPE_EDGE:
00210 printf("METIS_IPTYPE_EDGE\n");
00211 break;
00212 case METIS_IPTYPE_NODE:
00213 printf("METIS_IPTYPE_NODE\n");
00214 break;
00215 case METIS_IPTYPE_METISRB:
00216 printf("METIS_IPTYPE_METISRB\n");
00217 break;
00218 default:
00219 printf("Unknown!\n");
00220 }
00221
00222 printf(" Refinement type: ");
00223 switch (ctrl->rtype) {
00224 case METIS_RTYPE_FM:
00225 printf("METIS_RTYPE_FM\n");
00226 break;
00227 case METIS_RTYPE_GREEDY:
00228 printf("METIS_RTYPE_GREEDY\n");
00229 break;
00230 case METIS_RTYPE_SEP2SIDED:
00231 printf("METIS_RTYPE_SEP2SIDED\n");
00232 break;
00233 case METIS_RTYPE_SEP1SIDED:
00234 printf("METIS_RTYPE_SEP1SIDED\n");
00235 break;
00236 default:
00237 printf("Unknown!\n");
00238 }
00239
00240 printf(" Perform a 2-hop matching: %s\n", (ctrl->no2hop ? "Yes" : "No"));
00241
00242 printf(" Number of balancing constraints: %"PRIDX"\n", ctrl->ncon);
00243 printf(" Number of refinement iterations: %"PRIDX"\n", ctrl->niter);
00244 printf(" Random number seed: %"PRIDX"\n", ctrl->seed);
00245
00246 if (ctrl->optype == METIS_OP_OMETIS) {
00247 printf(" Number of separators: %"PRIDX"\n", ctrl->nseps);
00248 printf(" Compress graph prior to ordering: %s\n", (ctrl->compress ? "Yes" : "No"));
00249 printf(" Detect & order connected components separately: %s\n", (ctrl->ccorder ? "Yes" : "No"));
00250 printf(" Prunning factor for high degree vertices: %"PRREAL"\n", ctrl->pfactor);
00251 }
00252 else {
00253 printf(" Number of partitions: %"PRIDX"\n", ctrl->nparts);
00254 printf(" Number of cuts: %"PRIDX"\n", ctrl->ncuts);
00255 printf(" User-supplied ufactor: %"PRIDX"\n", ctrl->ufactor);
00256
00257 if (ctrl->optype == METIS_OP_KMETIS) {
00258 printf(" Minimize connectivity: %s\n", (ctrl->minconn ? "Yes" : "No"));
00259 printf(" Create contigous partitions: %s\n", (ctrl->contig ? "Yes" : "No"));
00260 }
00261
00262 modnum = (ctrl->ncon==1 ? 5 : (ctrl->ncon==2 ? 3 : (ctrl->ncon==3 ? 2 : 1)));
00263 printf(" Target partition weights: ");
00264 for (i=0; i<ctrl->nparts; i++) {
00265 if (i%modnum == 0)
00266 printf("\n ");
00267 printf("%4"PRIDX"=[", i);
00268 for (j=0; j<ctrl->ncon; j++)
00269 printf("%s%.2e", (j==0 ? "" : " "), (double)ctrl->tpwgts[i*ctrl->ncon+j]);
00270 printf("]");
00271 }
00272 printf("\n");
00273 }
00274
00275 printf(" Allowed maximum load imbalance: ");
00276 for (i=0; i<ctrl->ncon; i++)
00277 printf("%.3"PRREAL" ", ctrl->ubfactors[i]);
00278 printf("\n");
00279
00280 printf("\n");
00281 }
00282
00283
00284
00286
00287 int CheckParams(ctrl_t *ctrl)
00288 {
00289 idx_t i, j;
00290 real_t sum;
00291 mdbglvl_et dbglvl=METIS_DBG_INFO;
00292
00293 switch (ctrl->optype) {
00294 case METIS_OP_PMETIS:
00295 if (ctrl->objtype != METIS_OBJTYPE_CUT) {
00296 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect objective type.\n"));
00297 return 0;
00298 }
00299 if (ctrl->ctype != METIS_CTYPE_RM && ctrl->ctype != METIS_CTYPE_SHEM) {
00300 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect coarsening scheme.\n"));
00301 return 0;
00302 }
00303 if (ctrl->iptype != METIS_IPTYPE_GROW && ctrl->iptype != METIS_IPTYPE_RANDOM) {
00304 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect initial partitioning scheme.\n"));
00305 return 0;
00306 }
00307 if (ctrl->rtype != METIS_RTYPE_FM) {
00308 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect refinement scheme.\n"));
00309 return 0;
00310 }
00311 if (ctrl->ncuts <= 0) {
00312 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncuts.\n"));
00313 return 0;
00314 }
00315 if (ctrl->niter <= 0) {
00316 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect niter.\n"));
00317 return 0;
00318 }
00319 if (ctrl->ufactor <= 0) {
00320 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ufactor.\n"));
00321 return 0;
00322 }
00323 if (ctrl->numflag != 0 && ctrl->numflag != 1) {
00324 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect numflag.\n"));
00325 return 0;
00326 }
00327 if (ctrl->nparts <= 0) {
00328 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nparts.\n"));
00329 return 0;
00330 }
00331 if (ctrl->ncon <= 0) {
00332 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncon.\n"));
00333 return 0;
00334 }
00335
00336 for (i=0; i<ctrl->ncon; i++) {
00337 sum = rsum(ctrl->nparts, ctrl->tpwgts+i, ctrl->ncon);
00338 if (sum < 0.99 || sum > 1.01) {
00339 IFSET(dbglvl, METIS_DBG_INFO,
00340 printf("Input Error: Incorrect sum of %"PRREAL" for tpwgts for constraint %"PRIDX".\n", sum, i));
00341 return 0;
00342 }
00343 }
00344 for (i=0; i<ctrl->ncon; i++) {
00345 for (j=0; j<ctrl->nparts; j++) {
00346 if (ctrl->tpwgts[j*ctrl->ncon+i] <= 0.0) {
00347 IFSET(dbglvl, METIS_DBG_INFO,
00348 printf("Input Error: Incorrect tpwgts for partition %"PRIDX" and constraint %"PRIDX".\n", j, i));
00349 return 0;
00350 }
00351 }
00352 }
00353
00354 for (i=0; i<ctrl->ncon; i++) {
00355 if (ctrl->ubfactors[i] <= 1.0) {
00356 IFSET(dbglvl, METIS_DBG_INFO,
00357 printf("Input Error: Incorrect ubfactor for constraint %"PRIDX".\n", i));
00358 return 0;
00359 }
00360 }
00361
00362 break;
00363
00364 case METIS_OP_KMETIS:
00365 if (ctrl->objtype != METIS_OBJTYPE_CUT && ctrl->objtype != METIS_OBJTYPE_VOL) {
00366 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect objective type.\n"));
00367 return 0;
00368 }
00369 if (ctrl->ctype != METIS_CTYPE_RM && ctrl->ctype != METIS_CTYPE_SHEM) {
00370 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect coarsening scheme.\n"));
00371 return 0;
00372 }
00373 if (ctrl->iptype != METIS_IPTYPE_METISRB) {
00374 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect initial partitioning scheme.\n"));
00375 return 0;
00376 }
00377 if (ctrl->rtype != METIS_RTYPE_GREEDY) {
00378 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect refinement scheme.\n"));
00379 return 0;
00380 }
00381 if (ctrl->ncuts <= 0) {
00382 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncuts.\n"));
00383 return 0;
00384 }
00385 if (ctrl->niter <= 0) {
00386 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect niter.\n"));
00387 return 0;
00388 }
00389 if (ctrl->ufactor <= 0) {
00390 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ufactor.\n"));
00391 return 0;
00392 }
00393 if (ctrl->numflag != 0 && ctrl->numflag != 1) {
00394 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect numflag.\n"));
00395 return 0;
00396 }
00397 if (ctrl->nparts <= 0) {
00398 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nparts.\n"));
00399 return 0;
00400 }
00401 if (ctrl->ncon <= 0) {
00402 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncon.\n"));
00403 return 0;
00404 }
00405 if (ctrl->contig != 0 && ctrl->contig != 1) {
00406 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect contig.\n"));
00407 return 0;
00408 }
00409 if (ctrl->minconn != 0 && ctrl->minconn != 1) {
00410 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect minconn.\n"));
00411 return 0;
00412 }
00413
00414 for (i=0; i<ctrl->ncon; i++) {
00415 sum = rsum(ctrl->nparts, ctrl->tpwgts+i, ctrl->ncon);
00416 if (sum < 0.99 || sum > 1.01) {
00417 IFSET(dbglvl, METIS_DBG_INFO,
00418 printf("Input Error: Incorrect sum of %"PRREAL" for tpwgts for constraint %"PRIDX".\n", sum, i));
00419 return 0;
00420 }
00421 }
00422 for (i=0; i<ctrl->ncon; i++) {
00423 for (j=0; j<ctrl->nparts; j++) {
00424 if (ctrl->tpwgts[j*ctrl->ncon+i] <= 0.0) {
00425 IFSET(dbglvl, METIS_DBG_INFO,
00426 printf("Input Error: Incorrect tpwgts for partition %"PRIDX" and constraint %"PRIDX".\n", j, i));
00427 return 0;
00428 }
00429 }
00430 }
00431
00432 for (i=0; i<ctrl->ncon; i++) {
00433 if (ctrl->ubfactors[i] <= 1.0) {
00434 IFSET(dbglvl, METIS_DBG_INFO,
00435 printf("Input Error: Incorrect ubfactor for constraint %"PRIDX".\n", i));
00436 return 0;
00437 }
00438 }
00439
00440 break;
00441
00442
00443
00444 case METIS_OP_OMETIS:
00445 if (ctrl->objtype != METIS_OBJTYPE_NODE) {
00446 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect objective type.\n"));
00447 return 0;
00448 }
00449 if (ctrl->ctype != METIS_CTYPE_RM && ctrl->ctype != METIS_CTYPE_SHEM) {
00450 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect coarsening scheme.\n"));
00451 return 0;
00452 }
00453 if (ctrl->iptype != METIS_IPTYPE_EDGE && ctrl->iptype != METIS_IPTYPE_NODE) {
00454 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect initial partitioning scheme.\n"));
00455 return 0;
00456 }
00457 if (ctrl->rtype != METIS_RTYPE_SEP1SIDED && ctrl->rtype != METIS_RTYPE_SEP2SIDED) {
00458 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect refinement scheme.\n"));
00459 return 0;
00460 }
00461 if (ctrl->nseps <= 0) {
00462 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nseps.\n"));
00463 return 0;
00464 }
00465 if (ctrl->niter <= 0) {
00466 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect niter.\n"));
00467 return 0;
00468 }
00469 if (ctrl->ufactor <= 0) {
00470 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ufactor.\n"));
00471 return 0;
00472 }
00473 if (ctrl->numflag != 0 && ctrl->numflag != 1) {
00474 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect numflag.\n"));
00475 return 0;
00476 }
00477 if (ctrl->nparts != 3) {
00478 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nparts.\n"));
00479 return 0;
00480 }
00481 if (ctrl->ncon != 1) {
00482 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncon.\n"));
00483 return 0;
00484 }
00485 if (ctrl->compress != 0 && ctrl->compress != 1) {
00486 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect compress.\n"));
00487 return 0;
00488 }
00489 if (ctrl->ccorder != 0 && ctrl->ccorder != 1) {
00490 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ccorder.\n"));
00491 return 0;
00492 }
00493 if (ctrl->pfactor < 0.0 ) {
00494 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect pfactor.\n"));
00495 return 0;
00496 }
00497
00498 for (i=0; i<ctrl->ncon; i++) {
00499 if (ctrl->ubfactors[i] <= 1.0) {
00500 IFSET(dbglvl, METIS_DBG_INFO,
00501 printf("Input Error: Incorrect ubfactor for constraint %"PRIDX".\n", i));
00502 return 0;
00503 }
00504 }
00505
00506 break;
00507
00508 default:
00509 IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect optype\n"));
00510 return 0;
00511 }
00512
00513 return 1;
00514 }
00515
00516
00517
00519
00520 void FreeCtrl(ctrl_t **r_ctrl)
00521 {
00522 ctrl_t *ctrl = *r_ctrl;
00523
00524 FreeWorkSpace(ctrl);
00525
00526 gk_free((void **)&ctrl->tpwgts, &ctrl->pijbm,
00527 &ctrl->ubfactors, &ctrl->maxvwgt, &ctrl, LTERM);
00528
00529 *r_ctrl = NULL;
00530 }
00531
00532