options.c
Go to the documentation of this file.
1 
11 #include "metislib.h"
12 
13 
14 /*************************************************************************/
16 /*************************************************************************/
19 {
20  idx_t i, j;
21  ctrl_t *ctrl;
22 
23  ctrl = (ctrl_t *)gk_malloc(sizeof(ctrl_t), "SetupCtrl: ctrl");
24 
25  memset((void *)ctrl, 0, sizeof(ctrl_t));
26 
27  switch (optype) {
28  case METIS_OP_PMETIS:
29  ctrl->objtype = GETOPTION(options, METIS_OPTION_OBJTYPE, METIS_OBJTYPE_CUT);
30  ctrl->rtype = METIS_RTYPE_FM;
31  ctrl->ncuts = GETOPTION(options, METIS_OPTION_NCUTS, 1);
32  ctrl->niter = GETOPTION(options, METIS_OPTION_NITER, 10);
33 
34  if (ncon == 1) {
35  ctrl->iptype = GETOPTION(options, METIS_OPTION_IPTYPE, METIS_IPTYPE_GROW);
36  ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, PMETIS_DEFAULT_UFACTOR);
37  ctrl->CoarsenTo = 20;
38  }
39  else {
40  ctrl->iptype = GETOPTION(options, METIS_OPTION_IPTYPE, METIS_IPTYPE_RANDOM);
41  ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, MCPMETIS_DEFAULT_UFACTOR);
42  ctrl->CoarsenTo = 100;
43  }
44 
45  break;
46 
47 
48  case METIS_OP_KMETIS:
49  ctrl->objtype = GETOPTION(options, METIS_OPTION_OBJTYPE, METIS_OBJTYPE_CUT);
50  ctrl->iptype = METIS_IPTYPE_METISRB;
51  ctrl->rtype = METIS_RTYPE_GREEDY;
52  ctrl->ncuts = GETOPTION(options, METIS_OPTION_NCUTS, 1);
53  ctrl->niter = GETOPTION(options, METIS_OPTION_NITER, 10);
54  ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, KMETIS_DEFAULT_UFACTOR);
55  ctrl->minconn = GETOPTION(options, METIS_OPTION_MINCONN, 0);
56  ctrl->contig = GETOPTION(options, METIS_OPTION_CONTIG, 0);
57  break;
58 
59 
60  case METIS_OP_OMETIS:
61  ctrl->objtype = GETOPTION(options, METIS_OPTION_OBJTYPE, METIS_OBJTYPE_NODE);
62  ctrl->rtype = GETOPTION(options, METIS_OPTION_RTYPE, METIS_RTYPE_SEP1SIDED);
63  ctrl->iptype = GETOPTION(options, METIS_OPTION_IPTYPE, METIS_IPTYPE_EDGE);
64  ctrl->nseps = GETOPTION(options, METIS_OPTION_NSEPS, 1);
65  ctrl->niter = GETOPTION(options, METIS_OPTION_NITER, 10);
66  ctrl->ufactor = GETOPTION(options, METIS_OPTION_UFACTOR, OMETIS_DEFAULT_UFACTOR);
67  ctrl->compress = GETOPTION(options, METIS_OPTION_COMPRESS, 1);
68  ctrl->ccorder = GETOPTION(options, METIS_OPTION_CCORDER, 0);
69  ctrl->pfactor = 0.1*GETOPTION(options, METIS_OPTION_PFACTOR, 0);
70 
71  ctrl->CoarsenTo = 100;
72  break;
73 
74  default:
75  gk_errexit(SIGERR, "Unknown optype of %d\n", optype);
76  }
77 
78  /* common options */
79  ctrl->ctype = GETOPTION(options, METIS_OPTION_CTYPE, METIS_CTYPE_SHEM);
80  ctrl->no2hop = GETOPTION(options, METIS_OPTION_NO2HOP, 0);
81  ctrl->seed = GETOPTION(options, METIS_OPTION_SEED, -1);
82  ctrl->dbglvl = GETOPTION(options, METIS_OPTION_DBGLVL, 0);
83  ctrl->numflag = GETOPTION(options, METIS_OPTION_NUMBERING, 0);
84 
85  /* set non-option information */
86  ctrl->optype = optype;
87  ctrl->ncon = ncon;
88  ctrl->nparts = nparts;
89  ctrl->maxvwgt = ismalloc(ncon, 0, "SetupCtrl: maxvwgt");
90 
91  /* setup the target partition weights */
92  if (ctrl->optype != METIS_OP_OMETIS) {
93  ctrl->tpwgts = rmalloc(nparts*ncon, "SetupCtrl: ctrl->tpwgts");
94  if (tpwgts) {
95  rcopy(nparts*ncon, tpwgts, ctrl->tpwgts);
96  }
97  else {
98  for (i=0; i<nparts; i++) {
99  for (j=0; j<ncon; j++)
100  ctrl->tpwgts[i*ncon+j] = 1.0/nparts;
101  }
102  }
103  }
104  else { /* METIS_OP_OMETIS */
105  /* this is required to allow the pijbm to be defined properly for
106  the edge-based refinement during initial partitioning */
107  ctrl->tpwgts = rsmalloc(2, .5, "SetupCtrl: ctrl->tpwgts");
108  }
109 
110 
111  /* setup the ubfactors */
112  ctrl->ubfactors = rsmalloc(ctrl->ncon, I2RUBFACTOR(ctrl->ufactor), "SetupCtrl: ubfactors");
113  if (ubvec)
114  rcopy(ctrl->ncon, ubvec, ctrl->ubfactors);
115  for (i=0; i<ctrl->ncon; i++)
116  ctrl->ubfactors[i] += 0.0000499;
117 
118  /* Allocate memory for balance multipliers.
119  Note that for PMETIS/OMETIS routines the memory allocated is more
120  than required as balance multipliers for 2 parts is sufficient. */
121  ctrl->pijbm = rmalloc(nparts*ncon, "SetupCtrl: ctrl->pijbm");
122 
123  InitRandom(ctrl->seed);
124 
125  IFSET(ctrl->dbglvl, METIS_DBG_INFO, PrintCtrl(ctrl));
126 
127  if (!CheckParams(ctrl)) {
128  FreeCtrl(&ctrl);
129  return NULL;
130  }
131  else {
132  return ctrl;
133  }
134 }
135 
136 
137 /*************************************************************************/
139 /*************************************************************************/
141 {
142  idx_t i, j;
143 
144  for (i=0; i<ctrl->nparts; i++) {
145  for (j=0; j<graph->ncon; j++)
146  ctrl->pijbm[i*graph->ncon+j] = graph->invtvwgt[j]/ctrl->tpwgts[i*graph->ncon+j];
147  }
148 }
149 
150 
151 /*************************************************************************/
153 /*************************************************************************/
155 {
156  idx_t i, j;
157 
158  for (i=0; i<2; i++) {
159  for (j=0; j<graph->ncon; j++)
160  ctrl->pijbm[i*graph->ncon+j] = graph->invtvwgt[j]/tpwgts[i*graph->ncon+j];
161  }
162 }
163 
164 
165 /*************************************************************************/
167 /*************************************************************************/
168 void PrintCtrl(ctrl_t *ctrl)
169 {
170  idx_t i, j, modnum;
171 
172  printf(" Runtime parameters:\n");
173 
174  printf(" Objective type: ");
175  switch (ctrl->objtype) {
176  case METIS_OBJTYPE_CUT:
177  printf("METIS_OBJTYPE_CUT\n");
178  break;
179  case METIS_OBJTYPE_VOL:
180  printf("METIS_OBJTYPE_VOL\n");
181  break;
182  case METIS_OBJTYPE_NODE:
183  printf("METIS_OBJTYPE_NODE\n");
184  break;
185  default:
186  printf("Unknown!\n");
187  }
188 
189  printf(" Coarsening type: ");
190  switch (ctrl->ctype) {
191  case METIS_CTYPE_RM:
192  printf("METIS_CTYPE_RM\n");
193  break;
194  case METIS_CTYPE_SHEM:
195  printf("METIS_CTYPE_SHEM\n");
196  break;
197  default:
198  printf("Unknown!\n");
199  }
200 
201  printf(" Initial partitioning type: ");
202  switch (ctrl->iptype) {
203  case METIS_IPTYPE_GROW:
204  printf("METIS_IPTYPE_GROW\n");
205  break;
206  case METIS_IPTYPE_RANDOM:
207  printf("METIS_IPTYPE_RANDOM\n");
208  break;
209  case METIS_IPTYPE_EDGE:
210  printf("METIS_IPTYPE_EDGE\n");
211  break;
212  case METIS_IPTYPE_NODE:
213  printf("METIS_IPTYPE_NODE\n");
214  break;
216  printf("METIS_IPTYPE_METISRB\n");
217  break;
218  default:
219  printf("Unknown!\n");
220  }
221 
222  printf(" Refinement type: ");
223  switch (ctrl->rtype) {
224  case METIS_RTYPE_FM:
225  printf("METIS_RTYPE_FM\n");
226  break;
227  case METIS_RTYPE_GREEDY:
228  printf("METIS_RTYPE_GREEDY\n");
229  break;
231  printf("METIS_RTYPE_SEP2SIDED\n");
232  break;
234  printf("METIS_RTYPE_SEP1SIDED\n");
235  break;
236  default:
237  printf("Unknown!\n");
238  }
239 
240  printf(" Perform a 2-hop matching: %s\n", (ctrl->no2hop ? "Yes" : "No"));
241 
242  printf(" Number of balancing constraints: %"PRIDX"\n", ctrl->ncon);
243  printf(" Number of refinement iterations: %"PRIDX"\n", ctrl->niter);
244  printf(" Random number seed: %"PRIDX"\n", ctrl->seed);
245 
246  if (ctrl->optype == METIS_OP_OMETIS) {
247  printf(" Number of separators: %"PRIDX"\n", ctrl->nseps);
248  printf(" Compress graph prior to ordering: %s\n", (ctrl->compress ? "Yes" : "No"));
249  printf(" Detect & order connected components separately: %s\n", (ctrl->ccorder ? "Yes" : "No"));
250  printf(" Prunning factor for high degree vertices: %"PRREAL"\n", ctrl->pfactor);
251  }
252  else {
253  printf(" Number of partitions: %"PRIDX"\n", ctrl->nparts);
254  printf(" Number of cuts: %"PRIDX"\n", ctrl->ncuts);
255  printf(" User-supplied ufactor: %"PRIDX"\n", ctrl->ufactor);
256 
257  if (ctrl->optype == METIS_OP_KMETIS) {
258  printf(" Minimize connectivity: %s\n", (ctrl->minconn ? "Yes" : "No"));
259  printf(" Create contigous partitions: %s\n", (ctrl->contig ? "Yes" : "No"));
260  }
261 
262  modnum = (ctrl->ncon==1 ? 5 : (ctrl->ncon==2 ? 3 : (ctrl->ncon==3 ? 2 : 1)));
263  printf(" Target partition weights: ");
264  for (i=0; i<ctrl->nparts; i++) {
265  if (i%modnum == 0)
266  printf("\n ");
267  printf("%4"PRIDX"=[", i);
268  for (j=0; j<ctrl->ncon; j++)
269  printf("%s%.2e", (j==0 ? "" : " "), (double)ctrl->tpwgts[i*ctrl->ncon+j]);
270  printf("]");
271  }
272  printf("\n");
273  }
274 
275  printf(" Allowed maximum load imbalance: ");
276  for (i=0; i<ctrl->ncon; i++)
277  printf("%.3"PRREAL" ", ctrl->ubfactors[i]);
278  printf("\n");
279 
280  printf("\n");
281 }
282 
283 
284 /*************************************************************************/
286 /*************************************************************************/
287 int CheckParams(ctrl_t *ctrl)
288 {
289  idx_t i, j;
290  real_t sum;
291  mdbglvl_et dbglvl=METIS_DBG_INFO;
292 
293  switch (ctrl->optype) {
294  case METIS_OP_PMETIS:
295  if (ctrl->objtype != METIS_OBJTYPE_CUT) {
296  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect objective type.\n"));
297  return 0;
298  }
299  if (ctrl->ctype != METIS_CTYPE_RM && ctrl->ctype != METIS_CTYPE_SHEM) {
300  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect coarsening scheme.\n"));
301  return 0;
302  }
303  if (ctrl->iptype != METIS_IPTYPE_GROW && ctrl->iptype != METIS_IPTYPE_RANDOM) {
304  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect initial partitioning scheme.\n"));
305  return 0;
306  }
307  if (ctrl->rtype != METIS_RTYPE_FM) {
308  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect refinement scheme.\n"));
309  return 0;
310  }
311  if (ctrl->ncuts <= 0) {
312  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncuts.\n"));
313  return 0;
314  }
315  if (ctrl->niter <= 0) {
316  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect niter.\n"));
317  return 0;
318  }
319  if (ctrl->ufactor <= 0) {
320  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ufactor.\n"));
321  return 0;
322  }
323  if (ctrl->numflag != 0 && ctrl->numflag != 1) {
324  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect numflag.\n"));
325  return 0;
326  }
327  if (ctrl->nparts <= 0) {
328  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nparts.\n"));
329  return 0;
330  }
331  if (ctrl->ncon <= 0) {
332  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncon.\n"));
333  return 0;
334  }
335 
336  for (i=0; i<ctrl->ncon; i++) {
337  sum = rsum(ctrl->nparts, ctrl->tpwgts+i, ctrl->ncon);
338  if (sum < 0.99 || sum > 1.01) {
339  IFSET(dbglvl, METIS_DBG_INFO,
340  printf("Input Error: Incorrect sum of %"PRREAL" for tpwgts for constraint %"PRIDX".\n", sum, i));
341  return 0;
342  }
343  }
344  for (i=0; i<ctrl->ncon; i++) {
345  for (j=0; j<ctrl->nparts; j++) {
346  if (ctrl->tpwgts[j*ctrl->ncon+i] <= 0.0) {
347  IFSET(dbglvl, METIS_DBG_INFO,
348  printf("Input Error: Incorrect tpwgts for partition %"PRIDX" and constraint %"PRIDX".\n", j, i));
349  return 0;
350  }
351  }
352  }
353 
354  for (i=0; i<ctrl->ncon; i++) {
355  if (ctrl->ubfactors[i] <= 1.0) {
356  IFSET(dbglvl, METIS_DBG_INFO,
357  printf("Input Error: Incorrect ubfactor for constraint %"PRIDX".\n", i));
358  return 0;
359  }
360  }
361 
362  break;
363 
364  case METIS_OP_KMETIS:
365  if (ctrl->objtype != METIS_OBJTYPE_CUT && ctrl->objtype != METIS_OBJTYPE_VOL) {
366  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect objective type.\n"));
367  return 0;
368  }
369  if (ctrl->ctype != METIS_CTYPE_RM && ctrl->ctype != METIS_CTYPE_SHEM) {
370  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect coarsening scheme.\n"));
371  return 0;
372  }
373  if (ctrl->iptype != METIS_IPTYPE_METISRB) {
374  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect initial partitioning scheme.\n"));
375  return 0;
376  }
377  if (ctrl->rtype != METIS_RTYPE_GREEDY) {
378  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect refinement scheme.\n"));
379  return 0;
380  }
381  if (ctrl->ncuts <= 0) {
382  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncuts.\n"));
383  return 0;
384  }
385  if (ctrl->niter <= 0) {
386  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect niter.\n"));
387  return 0;
388  }
389  if (ctrl->ufactor <= 0) {
390  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ufactor.\n"));
391  return 0;
392  }
393  if (ctrl->numflag != 0 && ctrl->numflag != 1) {
394  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect numflag.\n"));
395  return 0;
396  }
397  if (ctrl->nparts <= 0) {
398  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nparts.\n"));
399  return 0;
400  }
401  if (ctrl->ncon <= 0) {
402  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncon.\n"));
403  return 0;
404  }
405  if (ctrl->contig != 0 && ctrl->contig != 1) {
406  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect contig.\n"));
407  return 0;
408  }
409  if (ctrl->minconn != 0 && ctrl->minconn != 1) {
410  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect minconn.\n"));
411  return 0;
412  }
413 
414  for (i=0; i<ctrl->ncon; i++) {
415  sum = rsum(ctrl->nparts, ctrl->tpwgts+i, ctrl->ncon);
416  if (sum < 0.99 || sum > 1.01) {
417  IFSET(dbglvl, METIS_DBG_INFO,
418  printf("Input Error: Incorrect sum of %"PRREAL" for tpwgts for constraint %"PRIDX".\n", sum, i));
419  return 0;
420  }
421  }
422  for (i=0; i<ctrl->ncon; i++) {
423  for (j=0; j<ctrl->nparts; j++) {
424  if (ctrl->tpwgts[j*ctrl->ncon+i] <= 0.0) {
425  IFSET(dbglvl, METIS_DBG_INFO,
426  printf("Input Error: Incorrect tpwgts for partition %"PRIDX" and constraint %"PRIDX".\n", j, i));
427  return 0;
428  }
429  }
430  }
431 
432  for (i=0; i<ctrl->ncon; i++) {
433  if (ctrl->ubfactors[i] <= 1.0) {
434  IFSET(dbglvl, METIS_DBG_INFO,
435  printf("Input Error: Incorrect ubfactor for constraint %"PRIDX".\n", i));
436  return 0;
437  }
438  }
439 
440  break;
441 
442 
443 
444  case METIS_OP_OMETIS:
445  if (ctrl->objtype != METIS_OBJTYPE_NODE) {
446  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect objective type.\n"));
447  return 0;
448  }
449  if (ctrl->ctype != METIS_CTYPE_RM && ctrl->ctype != METIS_CTYPE_SHEM) {
450  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect coarsening scheme.\n"));
451  return 0;
452  }
453  if (ctrl->iptype != METIS_IPTYPE_EDGE && ctrl->iptype != METIS_IPTYPE_NODE) {
454  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect initial partitioning scheme.\n"));
455  return 0;
456  }
457  if (ctrl->rtype != METIS_RTYPE_SEP1SIDED && ctrl->rtype != METIS_RTYPE_SEP2SIDED) {
458  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect refinement scheme.\n"));
459  return 0;
460  }
461  if (ctrl->nseps <= 0) {
462  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nseps.\n"));
463  return 0;
464  }
465  if (ctrl->niter <= 0) {
466  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect niter.\n"));
467  return 0;
468  }
469  if (ctrl->ufactor <= 0) {
470  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ufactor.\n"));
471  return 0;
472  }
473  if (ctrl->numflag != 0 && ctrl->numflag != 1) {
474  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect numflag.\n"));
475  return 0;
476  }
477  if (ctrl->nparts != 3) {
478  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect nparts.\n"));
479  return 0;
480  }
481  if (ctrl->ncon != 1) {
482  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ncon.\n"));
483  return 0;
484  }
485  if (ctrl->compress != 0 && ctrl->compress != 1) {
486  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect compress.\n"));
487  return 0;
488  }
489  if (ctrl->ccorder != 0 && ctrl->ccorder != 1) {
490  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect ccorder.\n"));
491  return 0;
492  }
493  if (ctrl->pfactor < 0.0 ) {
494  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect pfactor.\n"));
495  return 0;
496  }
497 
498  for (i=0; i<ctrl->ncon; i++) {
499  if (ctrl->ubfactors[i] <= 1.0) {
500  IFSET(dbglvl, METIS_DBG_INFO,
501  printf("Input Error: Incorrect ubfactor for constraint %"PRIDX".\n", i));
502  return 0;
503  }
504  }
505 
506  break;
507 
508  default:
509  IFSET(dbglvl, METIS_DBG_INFO, printf("Input Error: Incorrect optype\n"));
510  return 0;
511  }
512 
513  return 1;
514 }
515 
516 
517 /*************************************************************************/
519 /*************************************************************************/
520 void FreeCtrl(ctrl_t **r_ctrl)
521 {
522  ctrl_t *ctrl = *r_ctrl;
523 
524  FreeWorkSpace(ctrl);
525 
526  gk_free((void **)&ctrl->tpwgts, &ctrl->pijbm,
527  &ctrl->ubfactors, &ctrl->maxvwgt, &ctrl, LTERM);
528 
529  *r_ctrl = NULL;
530 }
531 
532 
#define I2RUBFACTOR(ufactor)
real_t * tpwgts
ctrl_t * SetupCtrl(moptype_et optype, idx_t *options, idx_t ncon, idx_t nparts, real_t *tpwgts, real_t *ubvec)
Definition: options.c:17
idx_t seed
mdbglvl_et
#define rcopy
Definition: gklib_rename.h:80
#define ismalloc
Definition: gklib_rename.h:68
idx_t ufactor
idx_t compress
idx_t idx_t idx_t idx_t idx_t idx_t idx_t * nparts
void gk_errexit(int signum, char *f_str,...)
Definition: error.c:77
idx_t ncuts
void SetupKWayBalMultipliers(ctrl_t *ctrl, graph_t *graph)
Definition: options.c:140
NonlinearFactorGraph graph
#define PRIDX
int CheckParams(ctrl_t *ctrl)
Definition: options.c:287
#define rsmalloc
Definition: gklib_rename.h:114
idx_t niter
idx_t ccorder
#define InitRandom
Definition: rename.h:246
#define rsum
Definition: gklib_rename.h:117
idx_t * maxvwgt
#define SIGERR
Definition: gk_defs.h:38
idx_t nseps
moptype_et optype
#define OMETIS_DEFAULT_UFACTOR
Definition: libmetis/defs.h:58
int32_t idx_t
void Setup2WayBalMultipliers(ctrl_t *ctrl, graph_t *graph, real_t *tpwgts)
Definition: options.c:154
real_t pfactor
idx_t nparts
idx_t no2hop
real_t * pijbm
idx_t ncon
idx_t idx_t idx_t idx_t idx_t idx_t idx_t real_t real_t * ubvec
float real_t
idx_t * ncon
#define NULL
Definition: ccolamd.c:609
idx_t idx_t idx_t idx_t idx_t idx_t idx_t real_t * tpwgts
mrtype_et rtype
idx_t ncon
void PrintCtrl(ctrl_t *ctrl)
Definition: options.c:168
real_t * ubfactors
void * gk_malloc(size_t nbytes, char *msg)
Definition: memory.c:140
void gk_free(void **ptr1,...)
Definition: memory.c:202
void FreeCtrl(ctrl_t **r_ctrl)
Definition: options.c:520
idx_t contig
#define FreeWorkSpace
Definition: rename.h:252
#define PRREAL
#define PMETIS_DEFAULT_UFACTOR
Definition: libmetis/defs.h:55
mctype_et ctype
#define IFSET(a, flag, cmd)
Definition: gk_macros.h:45
real_t * invtvwgt
moptype_et
idx_t minconn
#define MCPMETIS_DEFAULT_UFACTOR
Definition: libmetis/defs.h:56
std::ptrdiff_t j
mobjtype_et objtype
idx_t numflag
#define GETOPTION(options, idx, defval)
#define rmalloc
Definition: gklib_rename.h:93
#define LTERM
Definition: gk_defs.h:14
#define KMETIS_DEFAULT_UFACTOR
Definition: libmetis/defs.h:57
miptype_et iptype


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:34:59