Skip to main content

Documentation Index

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

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

Adaptive simulations change the mesh as computation proceeds — through local refinement, coarsening, or load changes — so that the initial partition becomes unbalanced over time. ParMETIS_V3_AdaptiveRepart repartitions an already-distributed graph while trying to minimize data movement between processes. Unlike a full repartition from scratch, it takes the existing partition into account and trades off partition quality against the cost of redistributing vertices.

When to repartition

Repartition when:
  • Adaptive mesh refinement or coarsening has shifted vertex counts significantly between processes.
  • Per-vertex computation costs (modeled by vwgt) have changed since the last partition, causing load imbalance to grow beyond the ubvec tolerance.
  • Communication volume has grown because the current partition no longer reflects the graph structure.
A common pattern is to check the maximum-to-average load ratio across processes after each time step and trigger repartitioning once it exceeds a threshold (e.g. 1.1 or 1.2).

The ipc2redist parameter

ipc2redist controls the tradeoff between inter-partition communication and redistribution cost.
ValueBehavior
Low (e.g. 1.0)Minimize data movement — keep vertices close to where they already are, even if partition quality suffers.
High (e.g. 1000.0)Maximize partition quality — redistribute aggressively to achieve the best cut, similar to a full repartition.
For most adaptive simulations a value between 10.0 and 100.0 provides a good balance. Start at 10.0 and increase if load imbalance remains high after repartitioning.

Coupling modes

The options[PMV3_OPTION_PSR] field controls whether the number of partitions must match the number of processes.
ConstantValueMeaning
PARMETIS_PSR_COUPLED1nparts == npes. Partition and process counts are locked together. Use this for the common case.
PARMETIS_PSR_UNCOUPLED2nparts may differ from npes. Allows over- or under-decomposition.

vsize: per-vertex communication volume

vsize is an optional array of length local_nvtxs that gives the communication volume associated with moving each vertex during redistribution. Pass NULL to treat all vertices as equal. Provide explicit values when vertices carry different amounts of data (e.g. variable-degree nodes with attached field data).

Workflow

1

Create an initial partition

Use ParMETIS_V3_PartKway (or load a checkpoint) to create the first partition. Record the part[] array on each process.
2

Run the simulation

Advance the simulation. Update the mesh and recompute vwgt to reflect the new per-vertex computational cost. Existing xadj and adjncy arrays should reflect the updated graph structure.
3

Update vertex weights

Set vwgt[i] to reflect the current compute cost of local vertex i. For a uniform-cost simulation you can keep all weights at 1. For adaptive simulations, heavier vertices might represent finer-resolution regions.
4

Call AdaptiveRepart

Pass the current part[] as input. After the call, part[] holds the updated partition and edgecut reports the new cut value.
options[0]               = 1;
options[PMV3_OPTION_PSR] = PARMETIS_PSR_COUPLED;

ParMETIS_V3_AdaptiveRepart(
    vtxdist, xadj, adjncy, vwgt,
    NULL,     /* vsize — NULL for equal vertex sizes */
    adjwgt,
    &wgtflag, &numflag, &ncon, &nparts,
    tpwgts, ubvec, &ipc2redist,
    options, &edgecut, part, &comm);
5

Redistribute data

Use the updated part[] to migrate vertex data to the new owning processes. Then continue the simulation.

Complete example

#include <stdio.h>
#include <stdlib.h>
#include <mpi.h>
#include <parmetis.h>

int main(int argc, char *argv[])
{
    MPI_Comm comm;
    int npes, mype;

    MPI_Init(&argc, &argv);
    MPI_Comm_dup(MPI_COMM_WORLD, &comm);
    MPI_Comm_size(comm, &npes);
    MPI_Comm_rank(comm, &mype);

    /*
     * Assume graph data (vtxdist, xadj, adjncy, vwgt, adjwgt)
     * has already been populated for the current mesh state.
     * Here we use a small ring graph for illustration.
     */
    idx_t vtxdist[3] = {0, 4, 8};   /* 8 vertices, 4 per process */
    idx_t nvtxs_local = vtxdist[mype+1] - vtxdist[mype];

    /* Minimal CSR ring graph (same as graph-partitioning example) */
    idx_t xadj_p0[]   = {0, 2, 4, 6, 8};
    idx_t adjncy_p0[] = {7, 1,  0, 2,  1, 3,  2, 4};
    idx_t xadj_p1[]   = {0, 2, 4, 6, 8};
    idx_t adjncy_p1[] = {3, 5,  4, 6,  5, 7,  6, 0};

    idx_t *xadj   = (mype == 0) ? xadj_p0   : xadj_p1;
    idx_t *adjncy = (mype == 0) ? adjncy_p0 : adjncy_p1;

    /* Vertex weights: simulate heavier load on process 0 after refinement */
    idx_t *vwgt   = (idx_t *)malloc(nvtxs_local * sizeof(idx_t));
    idx_t *adjwgt = (idx_t *)malloc(xadj[nvtxs_local] * sizeof(idx_t));
    for (idx_t i = 0; i < nvtxs_local; i++)
        vwgt[i] = (mype == 0) ? 3 : 1;    /* process 0 is 3x heavier */
    for (idx_t i = 0; i < xadj[nvtxs_local]; i++)
        adjwgt[i] = 1;

    idx_t  wgtflag = 3, numflag = 0, ncon = 1, nparts = 2;
    real_t tpwgts[2]  = {0.5f, 0.5f};
    real_t ubvec[1]   = {1.05f};
    real_t ipc2redist = 10.0f;   /* moderate redistribution tolerance */

    idx_t options[4] = {0};
    options[0]               = 1;
    options[PMV3_OPTION_PSR] = PARMETIS_PSR_COUPLED;

    /* Start from current partition (each vertex on its owning process) */
    idx_t *part = (idx_t *)malloc(nvtxs_local * sizeof(idx_t));
    for (idx_t i = 0; i < nvtxs_local; i++)
        part[i] = mype;

    idx_t edgecut;

    ParMETIS_V3_AdaptiveRepart(
        vtxdist, xadj, adjncy, vwgt,
        NULL,     /* vsize */
        adjwgt,
        &wgtflag, &numflag, &ncon, &nparts,
        tpwgts, ubvec, &ipc2redist,
        options, &edgecut, part, &comm);

    if (mype == 0)
        printf("After adaptive repartition, edge cut: %d\n", (int)edgecut);

    free(vwgt); free(adjwgt); free(part);
    MPI_Comm_free(&comm);
    MPI_Finalize();
    return 0;
}

Function signature

int ParMETIS_V3_AdaptiveRepart(
    idx_t *vtxdist, idx_t *xadj, idx_t *adjncy,
    idx_t *vwgt, idx_t *vsize, idx_t *adjwgt,
    idx_t *wgtflag, idx_t *numflag, idx_t *ncon, idx_t *nparts,
    real_t *tpwgts, real_t *ubvec, real_t *ipc2redist,
    idx_t *options, idx_t *edgecut, idx_t *part, MPI_Comm *comm);

Compilation

mpicc -o adaptive adaptive.c -lparmetis -lmetis
AdaptiveRepart reads the current contents of part[] as the starting partition. Passing an uninitialized or all-zero part[] is valid but may produce a suboptimal result because the algorithm has no useful migration history to minimize against.
If ipc2redist is set very high (e.g. 1000.0) the result will be close to calling PartKway from scratch. Use high values only when load imbalance is severe and data movement cost is acceptable.

Build docs developers (and LLMs) love