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.

METIS partitions finite element meshes by first converting them to a graph, partitioning the graph, and then deriving the other side’s partition. Two strategies are available: the dual approach treats elements as graph vertices (two elements are adjacent if they share at least ncommon nodes), and the nodal approach treats mesh nodes as graph vertices (two nodes are adjacent if they belong to the same element). The conversion functions METIS_MeshToDual and METIS_MeshToNodal can also be called directly when you need the graph representation itself.

Element-Node CSR Format

All mesh functions represent meshes using two arrays in a format analogous to graph CSR:
  • eptr — integer array of size ne + 1. Element i uses nodes eind[eptr[i] .. eptr[i+1]-1].
  • eind — flat array of node indices for all elements. Total length equals eptr[ne], which is the sum of the number of nodes over all elements.
Example: 2 triangular elements sharing an edge
  Element 0: nodes 0, 1, 2
  Element 1: nodes 1, 3, 2

eptr = [0, 3, 6]
eind = [0, 1, 2,   1, 3, 2]
Elements can have different numbers of nodes (mixed meshes), as long as eptr correctly marks each element’s slice of eind.

METIS_PartMeshDual

Partitions a mesh by constructing the dual graph (elements as vertices, shared faces as edges), partitioning it using METIS_PartGraphKway or METIS_PartGraphRecursive, and then inducing a node partition from the element partition.
int METIS_PartMeshDual(
    idx_t  *ne,
    idx_t  *nn,
    idx_t  *eptr,
    idx_t  *eind,
    idx_t  *vwgt,
    idx_t  *vsize,
    idx_t  *ncommon,
    idx_t  *nparts,
    real_t *tpwgts,
    idx_t  *options,
    idx_t  *objval,
    idx_t  *epart,
    idx_t  *npart
);
ne
idx_t *
required
Number of elements in the mesh.
nn
idx_t *
required
Number of nodes in the mesh.
eptr
idx_t *
required
Element-node CSR row-pointer array of size ne + 1.
eind
idx_t *
required
Element-node CSR column-index array. Contains the node indices of every element.
vwgt
idx_t *
Element weights array of size ne. Pass NULL for unit element weights.
vsize
idx_t *
Element sizes for communication-volume minimization, array of size ne. Pass NULL for unit sizes.
ncommon
idx_t *
required
Minimum number of nodes that two elements must share to be connected in the dual graph. For triangular/tetrahedral meshes with shared edges/faces, typical values are 2 (shared edge) or 3 (shared triangular face). Setting ncommon = 1 connects all elements that share any node.
nparts
idx_t *
required
Number of desired partitions.
tpwgts
real_t *
Target partition weights. Array of size nparts. tpwgts[i] is the desired fraction of the total element weight in partition i. Values must sum to 1.0. Pass NULL for equal partitioning.
options
idx_t *
required
Options array of size METIS_NOPTIONS initialized with METIS_SetDefaultOptions. The METIS_OPTION_PTYPE option selects between k-way (METIS_PTYPE_KWAY, default) and recursive bisection (METIS_PTYPE_RB).
objval
idx_t *
required
Output: objective value of the computed partition (edge cut or communication volume of the dual graph).
epart
idx_t *
required
Output: element partition array of size ne. epart[i] is the partition number assigned to element i. Must be pre-allocated by the caller.
npart
idx_t *
required
Output: node partition array of size nn. npart[i] is the partition number assigned to node i. Derived from the element partition by assigning each node to the partition that dominates its element neighborhood. Must be pre-allocated by the caller.

METIS_PartMeshNodal

Partitions a mesh by constructing the nodal graph (nodes as vertices, edges between nodes sharing an element), partitioning it, and inducing an element partition from the node partition.
int METIS_PartMeshNodal(
    idx_t  *ne,
    idx_t  *nn,
    idx_t  *eptr,
    idx_t  *eind,
    idx_t  *vwgt,
    idx_t  *vsize,
    idx_t  *nparts,
    real_t *tpwgts,
    idx_t  *options,
    idx_t  *objval,
    idx_t  *epart,
    idx_t  *npart
);
ne
idx_t *
required
Number of elements in the mesh.
nn
idx_t *
required
Number of nodes in the mesh.
eptr
idx_t *
required
Element-node CSR row-pointer array of size ne + 1.
eind
idx_t *
required
Element-node CSR column-index array.
vwgt
idx_t *
Node weights array of size nn. Pass NULL for unit node weights.
vsize
idx_t *
Node sizes for communication-volume minimization, array of size nn. Pass NULL for unit sizes.
nparts
idx_t *
required
Number of desired partitions.
tpwgts
real_t *
Target partition weights. Array of size nparts. Pass NULL for equal partitioning.
options
idx_t *
required
Options array initialized with METIS_SetDefaultOptions.
objval
idx_t *
required
Output: objective value of the computed partition.
epart
idx_t *
required
Output: element partition array of size ne. Derived from the node partition. Must be pre-allocated.
npart
idx_t *
required
Output: node partition array of size nn. Primary partitioning output. Must be pre-allocated.

METIS_MeshToDual

Converts a mesh to its dual graph. Returns the graph’s CSR arrays via output pointer-to-pointer parameters. The caller is responsible for freeing the returned arrays using METIS_Free.
int METIS_MeshToDual(
    idx_t  *ne,
    idx_t  *nn,
    idx_t  *eptr,
    idx_t  *eind,
    idx_t  *ncommon,
    idx_t  *numflag,
    idx_t **r_xadj,
    idx_t **r_adjncy
);
ne
idx_t *
required
Number of elements.
nn
idx_t *
required
Number of nodes.
eptr
idx_t *
required
Element-node CSR row-pointer array.
eind
idx_t *
required
Element-node CSR column-index array.
ncommon
idx_t *
required
Minimum number of shared nodes required for an edge in the dual graph.
numflag
idx_t *
required
Numbering flag. Pass 0 for 0-based numbering, 1 for 1-based (Fortran-style). The returned graph uses the same numbering.
r_xadj
idx_t **
required
Output: pointer to the allocated xadj array for the dual graph. METIS allocates this array. Free with METIS_Free(*r_xadj).
r_adjncy
idx_t **
required
Output: pointer to the allocated adjncy array for the dual graph. METIS allocates this array. Free with METIS_Free(*r_adjncy).

METIS_MeshToNodal

Converts a mesh to its nodal graph. Like METIS_MeshToDual, the returned arrays are allocated by METIS and must be freed by the caller with METIS_Free.
int METIS_MeshToNodal(
    idx_t  *ne,
    idx_t  *nn,
    idx_t  *eptr,
    idx_t  *eind,
    idx_t  *numflag,
    idx_t **r_xadj,
    idx_t **r_adjncy
);
ne
idx_t *
required
Number of elements.
nn
idx_t *
required
Number of nodes.
eptr
idx_t *
required
Element-node CSR row-pointer array.
eind
idx_t *
required
Element-node CSR column-index array.
numflag
idx_t *
required
Numbering flag: 0 for 0-based, 1 for 1-based.
r_xadj
idx_t **
required
Output: pointer to the allocated xadj array for the nodal graph. Free with METIS_Free(*r_xadj).
r_adjncy
idx_t **
required
Output: pointer to the allocated adjncy array for the nodal graph. Free with METIS_Free(*r_adjncy).
METIS_MeshToDual and METIS_MeshToNodal allocate the output arrays internally using the system malloc. You must free them by calling METIS_Free(ptr), not the standard free() directly (though they are equivalent in current versions). Do not use gk_free or any other METIS-internal allocator.

Complete Example: Partitioning a Triangular Mesh

#include <stdio.h>
#include <metis.h>

int main(void) {
    /*
     * Simple mesh: 4 triangles, 6 nodes
     *
     *   0 --- 1 --- 2
     *   |\ t0 |\ t2 |
     *   | \   | \   |
     *   |  \  |  \  |
     *   | t1\ | t3\  |
     *   3 --- 4 --- 5
     *
     * t0: nodes 0,1,4   t1: nodes 0,4,3
     * t2: nodes 1,2,5   t3: nodes 1,5,4
     */
    idx_t ne = 4, nn = 6;
    idx_t nparts = 2;
    idx_t ncommon = 2;   /* share an edge (2 nodes) to be adjacent */

    idx_t eptr[] = {0, 3, 6, 9, 12};
    idx_t eind[] = {
        0, 1, 4,    /* triangle 0 */
        0, 4, 3,    /* triangle 1 */
        1, 2, 5,    /* triangle 2 */
        1, 5, 4     /* triangle 3 */
    };

    idx_t epart[4], npart[6];
    idx_t objval;

    idx_t options[METIS_NOPTIONS];
    METIS_SetDefaultOptions(options);

    int ret = METIS_PartMeshDual(
        &ne, &nn,
        eptr, eind,
        NULL,      /* unit element weights */
        NULL,      /* no vsize */
        &ncommon,
        &nparts,
        NULL,      /* equal target weights */
        options,
        &objval,
        epart, npart
    );

    if (ret != METIS_OK) {
        fprintf(stderr, "METIS error %d\n", ret);
        return 1;
    }

    printf("Objective: %d\n", (int)objval);
    for (int i = 0; i < ne; i++)
        printf("  element %d -> part %d\n", i, (int)epart[i]);
    for (int i = 0; i < nn; i++)
        printf("  node %d    -> part %d\n", i, (int)npart[i]);

    return 0;
}

Using MeshToDual Directly

When you need the graph itself (for example, to inspect connectivity or to pass it to METIS_PartGraphKway with custom options):
idx_t *xadj = NULL, *adjncy = NULL;
idx_t numflag = 0;

int ret = METIS_MeshToDual(&ne, &nn, eptr, eind, &ncommon, &numflag,
                           &xadj, &adjncy);
if (ret != METIS_OK) { /* handle error */ }

/* Use xadj / adjncy ... */

/* Free when done — these were allocated by METIS */
METIS_Free(xadj);
METIS_Free(adjncy);

Dual vs Nodal: When to Use Each

Dual (PartMeshDual)Nodal (PartMeshNodal)
Primary partitionElementsNodes
Graph built fromElements sharing ncommon nodesNodes sharing an element
Typical useFEM solvers that work element-by-elementFEM solvers that assemble by node
ncommon parameterRequired (controls connectivity)Not applicable
Induced partitionNodes derived from elementsElements derived from nodes

Build docs developers (and LLMs) love