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.

Finite element applications distribute meshes across MPI processes as collections of elements, not vertices. ParMETIS provides two functions for this case: ParMETIS_V3_PartMeshKway partitions the elements directly, while ParMETIS_V3_Mesh2Dual converts a distributed element mesh to a dual graph so that ParMETIS_V3_PartKway can be used. Both approaches require describing the mesh with elmdist, eptr, and eind arrays, which are the element-level analogues of vtxdist, xadj, and adjncy.

Describing the distributed mesh

elmdist

elmdist is a replicated array of length npes + 1. elmdist[i] is the global index of the first element owned by process i. All processes must hold an identical copy.
/* 9 elements distributed across 3 processes: 3 elements each */
idx_t elmdist[4] = {0, 3, 6, 9};

eptr and eind

eptr and eind describe the node lists of local elements in CSR format. eptr[i+1] - eptr[i] is the number of nodes in local element i, and eind[eptr[i]..eptr[i+1]-1] holds their global node indices.
/* Process 0 owns 3 triangular elements, each with 3 nodes */
idx_t eptr[4]    = {0, 3, 6, 9};
idx_t eind[9]    = {0, 1, 3,   /* element 0: nodes 0, 1, 3 */
                    1, 4, 3,   /* element 1: nodes 1, 4, 3 */
                    1, 2, 4};  /* element 2: nodes 1, 2, 4 */

ncommonnodes

Two elements are considered adjacent when they share at least ncommonnodes nodes.
Element typeRecommended valueAdjacency kind
Triangles, quads2Edge adjacency
Tetrahedra3Face adjacency
Hexahedra4Face adjacency

Two approaches: PartMeshKway vs Mesh2Dual + PartKway

PartMeshKway is the simpler path. Call it directly when you only need a partition array for the elements. Mesh2Dual + PartKway gives more control. Mesh2Dual constructs the dual graph (one graph node per mesh element, edges between adjacent elements) and returns xadj/adjncy pointers that you then pass to PartKway. Use this when you need the dual graph for other purposes, or when you want to apply multi-constraint partitioning that PartMeshKway does not expose.
The xadj and adjncy arrays returned by Mesh2Dual are allocated by ParMETIS. You are responsible for freeing them with free() after use.

Complete example: PartMeshKway

#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);

    /*
     * 6 triangular elements distributed across 2 processes:
     *   process 0 owns elements 0-2
     *   process 1 owns elements 3-5
     *
     * Node layout (global indices):
     *   0 - 1 - 2
     *   |\ | /|
     *   3 - 4 - 5
     *   |\ | /|
     *   6 - 7 - 8
     */
    idx_t elmdist[3] = {0, 3, 6};

    /* Process 0: elements 0, 1, 2 (triangles) */
    idx_t eptr_p0[4]   = {0, 3, 6, 9};
    idx_t eind_p0[9]   = {0, 1, 3,
                          1, 4, 3,
                          1, 2, 4};

    /* Process 1: elements 3, 4, 5 (triangles) */
    idx_t eptr_p1[4]   = {0, 3, 6, 9};
    idx_t eind_p1[9]   = {3, 4, 6,
                          4, 7, 6,
                          4, 5, 7};

    idx_t *eptr = (mype == 0) ? eptr_p0 : eptr_p1;
    idx_t *eind = (mype == 0) ? eind_p0 : eind_p1;

    idx_t local_nelms = elmdist[mype+1] - elmdist[mype];

    /* No element weights */
    idx_t  wgtflag = 0;
    idx_t  numflag = 0;
    idx_t  ncon    = 1;
    idx_t  ncommonnodes = 2;   /* edge adjacency for triangles */
    idx_t  nparts  = 2;

    real_t tpwgts[2] = {0.5f, 0.5f};
    real_t ubvec[1]  = {1.05f};
    idx_t  options[3] = {0, 0, 0};

    idx_t  edgecut;
    idx_t *part = (idx_t *)malloc(local_nelms * sizeof(idx_t));

    ParMETIS_V3_PartMeshKway(
        elmdist, eptr, eind,
        NULL,           /* elmwgt — NULL because wgtflag=0 */
        &wgtflag, &numflag, &ncon, &ncommonnodes, &nparts,
        tpwgts, ubvec, options, &edgecut, part, &comm);

    if (mype == 0)
        printf("Edge cut: %d\n", (int)edgecut);

    for (int p = 0; p < npes; p++) {
        if (mype == p) {
            for (idx_t i = 0; i < local_nelms; i++)
                printf("process %d: element %d -> partition %d\n",
                       mype, (int)(elmdist[mype] + i), (int)part[i]);
        }
        MPI_Barrier(comm);
    }

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

Using Mesh2Dual

idx_t  numflag = 0, ncommonnodes = 2;
idx_t *xadj   = NULL;
idx_t *adjncy = NULL;

ParMETIS_V3_Mesh2Dual(
    elmdist, eptr, eind,
    &numflag, &ncommonnodes,
    &xadj, &adjncy, &comm);

/* Pass xadj and adjncy to PartKway ... */

/* Caller must free the returned arrays */
free(xadj);
free(adjncy);

Function signatures

int ParMETIS_V3_Mesh2Dual(
    idx_t *elmdist, idx_t *eptr, idx_t *eind,
    idx_t *numflag, idx_t *ncommonnodes,
    idx_t **xadj, idx_t **adjncy, MPI_Comm *comm);

int ParMETIS_V3_PartMeshKway(
    idx_t *elmdist, idx_t *eptr, idx_t *eind,
    idx_t *elmwgt,
    idx_t *wgtflag, idx_t *numflag, idx_t *ncon,
    idx_t *ncommonnodes, idx_t *nparts,
    real_t *tpwgts, real_t *ubvec, idx_t *options,
    idx_t *edgecut, idx_t *part, MPI_Comm *comm);

Compilation

mpicc -o meshpart meshpart.c -lparmetis -lmetis

Build docs developers (and LLMs) love