Skip to main content

Documentation Index

Fetch the complete documentation index at: https://mintlify.com/KarypisLab/METIS/llms.txt

Use this file to discover all available pages before exploring further.

Standard graph partitioning balances a single quantity—typically computational work—across partitions. In practice, parallel applications often involve multiple resources that must be balanced simultaneously: computation time, memory usage, communication bandwidth, or I/O load. METIS supports this through multi-constraint partitioning, where each vertex carries a vector of weights and the partitioner must balance every component of that weight vector across all partitions at the same time.

Enabling multi-constraint partitioning

Multi-constraint partitioning is enabled by setting ncon > 1 in the METIS_PartGraphKway or METIS_PartGraphRecursive call. The vwgt array then stores ncon weights per vertex, laid out as a flat interleaved array:
/* vertex i, constraint j → vwgt[i * ncon + j] */
idx_t vwgt[nvtxs * ncon];
For example, with nvtxs = 4 and ncon = 3:
vwgt index:   0  1  2    3  4  5    6  7  8    9 10 11
vertex:       ----0----  ----1----  ----2----  ----3----
constraint:    0  1  2    0  1  2    0  1  2    0  1  2
The number of constraints affects coarsening (the BetterVBalance function in mcutil.c prefers matchings that improve multi-dimensional balance), initial partitioning (multi-constraint bisection uses McRandomBisection and McGrowBisection), and refinement (balance checks test all ncon constraints simultaneously).
The sample file graphs/test.mgraph uses fmt=010 with ncon=2, giving each of its 766 vertices two weight values. This is a convenient file for testing multi-constraint behavior with gpmetis graphs/test.mgraph 4.

Target partition weights: tpwgts

The tpwgts array specifies the desired fraction of total weight that each partition should receive, independently per constraint. It is a flat array of size nparts * ncon:
/* partition i, constraint j → tpwgts[i * ncon + j] */
real_t tpwgts[nparts * ncon];
For each constraint j, the weights must sum to 1.0:
∑_i  tpwgts[i * ncon + j]  ==  1.0   for all j
Passing NULL for tpwgts is equivalent to setting every entry to 1.0 / nparts (equal partition weights for all constraints).

Example: 3 constraints, 4 partitions

Suppose you have 4 partitions and 3 constraints (compute, memory, I/O), and you want equal distribution across all:
idx_t nparts = 4;
idx_t ncon   = 3;

real_t tpwgts[4 * 3] = {
  /* part 0 */  0.25, 0.25, 0.25,
  /* part 1 */  0.25, 0.25, 0.25,
  /* part 2 */  0.25, 0.25, 0.25,
  /* part 3 */  0.25, 0.25, 0.25,
};
For a heterogeneous machine where partition 0 is twice as powerful as the others for compute (constraint 0) but equal for memory and I/O:
real_t tpwgts[4 * 3] = {
  /* part 0: 40% compute, 25% memory, 25% I/O */
  0.40, 0.25, 0.25,
  /* part 1: 20% compute, 25% memory, 25% I/O */
  0.20, 0.25, 0.25,
  /* part 2: 20% compute, 25% memory, 25% I/O */
  0.20, 0.25, 0.25,
  /* part 3: 20% compute, 25% memory, 25% I/O */
  0.20, 0.25, 0.25,
};
/* Each column (constraint) sums to 1.0 */

Imbalance tolerance: ubvec

The ubvec array specifies the maximum allowed load imbalance for each constraint independently. It has ncon elements:
real_t ubvec[ncon];
ubvec[j] is the multiplicative imbalance tolerance for constraint j. A value of 1.05 means each partition is allowed to hold up to 5% more than its target weight tpwgts[i*ncon+j] of constraint j’s total weight. The values must all be greater than 1.0. Passing NULL for ubvec sets the default tolerance: 1.001 when ncon == 1, and 1.01 when ncon > 1. Looser tolerances (larger ubvec values) allow METIS more freedom to improve cut quality at the cost of balance. Tighter tolerances force near-perfect balance but may produce higher edge cuts.

Complete C example with ncon=2

The following example uses the two-constraint structure of test.mgraph—each vertex has two weights—and partitions into 4 equal parts:
#include "metis.h"
#include <stdlib.h>

int main(void)
{
    idx_t nvtxs  = 766;
    idx_t ncon   = 2;   /* two vertex-weight constraints */
    idx_t nparts = 4;

    /* xadj and adjncy: standard CSR adjacency (built from test.mgraph) */
    idx_t *xadj   = /* ... */;
    idx_t *adjncy = /* ... */;

    /* vwgt: two weights per vertex, interleaved */
    idx_t *vwgt = /* ... */;  /* size: nvtxs * ncon */

    /* Equal target weights: each partition gets 1/4 of each constraint */
    real_t *tpwgts = NULL;   /* NULL → uniform 1/nparts */

    /* Allow 5% imbalance per constraint */
    real_t ubvec[2] = {1.05, 1.05};

    /* Initialize options to defaults */
    idx_t options[METIS_NOPTIONS];
    METIS_SetDefaultOptions(options);
    options[METIS_OPTION_OBJTYPE] = METIS_OBJTYPE_CUT;

    idx_t  objval;
    idx_t *part = malloc(nvtxs * sizeof(idx_t));

    int ret = METIS_PartGraphKway(
        &nvtxs, &ncon,
        xadj, adjncy,
        vwgt,   /* vertex weights: ncon values per vertex */
        NULL,   /* vsize: NULL → uniform sizes */
        NULL,   /* adjwgt: NULL → uniform edge weights */
        &nparts,
        tpwgts,
        ubvec,
        options,
        &objval,
        part
    );

    if (ret == METIS_OK) {
        /* part[i] in [0, nparts-1] for each vertex i */
    }

    free(part);
    return ret == METIS_OK ? 0 : 1;
}

When to use multi-constraint partitioning

Multi-constraint partitioning is useful whenever a single weight per vertex cannot adequately model the true computational requirements. Common scenarios:
Different processors have different compute speeds, memory capacities, or network bandwidths. Assign one constraint per resource type and set tpwgts to reflect the relative capabilities of each processor. Partition 0 might receive a larger share of compute weight if it is faster, while all partitions receive equal memory weight.
A mesh vertex may represent a point where multiple physical fields are computed (fluid, temperature, chemical species). Each field has different computational cost. Assigning one weight per field allows METIS to balance each field’s total work independently.
In GPU or accelerator programming, balancing computation alone may cause memory overflow on some devices. Adding a memory-usage constraint (e.g., bytes of working set per vertex) ensures that no partition exceeds device memory capacity.
A sparse matrix row may have different computational costs depending on whether it involves floating-point, integer, or boolean operations. Assigning one weight per operation type prevents load imbalance from mixed-type computations.

The UFACTOR option and balance interaction

METIS_OPTION_UFACTOR sets a global imbalance tolerance as an integer, used internally by METIS when ubvec is NULL. The relationship is:
ubvec[j] = 1.0 + UFACTOR / 1000.0
The default UFACTOR is 1 for METIS_PTYPE_RB (giving ubvec = 1.001) and 30 for METIS_PTYPE_KWAY (giving ubvec = 1.03). When you supply ubvec explicitly, UFACTOR is ignored for that call. For multi-constraint problems it is generally better to supply ubvec directly so you can control the tolerance per constraint, since different resource types may have very different sensitivity to imbalance.
Setting ubvec[j] too close to 1.0 (e.g., 1.001 with many constraints) can make it impossible for METIS to find a balanced partition, especially when the vertex weights are large integers. If METIS returns a partition that violates balance, try relaxing ubvec values slightly (e.g., 1.05 to 1.10) before debugging the input data.

Build docs developers (and LLMs) love