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.

ParMETIS_V3_PartMeshKway partitions a distributed finite element mesh without requiring an explicit dual graph construction step. Internally it builds the dual graph and applies the multilevel k-way algorithm in a single call, making it more convenient than calling ParMETIS_V3_Mesh2Dual followed by ParMETIS_V3_PartKway. Element adjacency is determined by the same ncommonnodes criterion used in Mesh2Dual. The resulting part array assigns each local element to a partition.
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);

Parameters

elmdist
idx_t*
Distribution of elements across processes. elmdist[i] is the global index of the first element on process i. Size npes + 1; must be identical on all processes.
eptr
idx_t*
Element connectivity row pointers. eptr[i] is the index into eind where the node list of local element i begins. Size local_nelems + 1.
eind
idx_t*
Element connectivity node indices. Contains the global mesh node indices for each local element. Size eptr[local_nelems].
elmwgt
idx_t*
Element weights. ncon values per element in interleaved order, analogous to vwgt for graphs. May be NULL when bit 1 of wgtflag is 0.
wgtflag
idx_t*
Bitmask controlling which weights are used. 0 = no weights, 2 = element weights only. Edge weights on the dual graph are not supported; bit 0 should be 0.
numflag
idx_t*
Indexing convention. 0 for C-style 0-based indexing; 1 for Fortran-style 1-based indexing.
ncon
idx_t*
Number of balancing constraints per element.
ncommonnodes
idx_t*
Minimum number of shared mesh nodes required for two elements to be considered adjacent. Use 2 for edge adjacency in 2D and 3 for face adjacency in 3D.
nparts
idx_t*
Desired number of partitions.
tpwgts
real_t*
Target partition weights. Array of nparts * ncon values. Values for each constraint must sum to 1.0. Pass NULL for equal-weight partitions.
ubvec
real_t*
Per-constraint imbalance tolerance. Array of ncon values greater than 1.0.
options
idx_t*
Algorithm options. If options[0] = 0, all defaults are used. If options[0] = 1, options[1] sets the debug level and options[2] sets the random seed.
edgecut
idx_t*
Output. Total weight of edges in the dual graph that cross partition boundaries.
part
idx_t*
Output. Partition assignment for each local element. part[i] is the partition ID of local element i. Caller must allocate at least local_nelems elements.
comm
MPI_Comm*
Pointer to the MPI communicator covering all participating processes.

Return value

Returns METIS_OK (1) on success. Any other value indicates an error.

Usage example

#include <parmetis.h>

int main(int argc, char *argv[]) {
    MPI_Init(&argc, &argv);
    MPI_Comm comm = MPI_COMM_WORLD;

    int npes, rank;
    MPI_Comm_size(comm, &npes);
    MPI_Comm_rank(comm, &rank);

    idx_t local_nelems = 8;
    idx_t elmdist[npes + 1];
    for (int i = 0; i <= npes; i++) elmdist[i] = i * local_nelems;

    /* Triangular mesh: 3 nodes per element */
    idx_t eptr[local_nelems + 1];
    idx_t eind[local_nelems * 3];
    for (int i = 0; i <= local_nelems; i++) eptr[i] = i * 3;
    /* populate eind with global node indices... */

    idx_t wgtflag      = 0;
    idx_t numflag      = 0;
    idx_t ncon         = 1;
    idx_t ncommonnodes = 2; /* edge adjacency for 2D */
    idx_t nparts       = npes;
    real_t tpwgts[nparts];
    for (int i = 0; i < nparts; i++) tpwgts[i] = 1.0f / nparts;
    real_t ubvec       = 1.05f;
    idx_t options[3]   = {0, 0, 0};
    idx_t edgecut;
    idx_t part[local_nelems];

    int ret = ParMETIS_V3_PartMeshKway(
        elmdist, eptr, eind,
        NULL,                    /* no element weights */
        &wgtflag, &numflag, &ncon, &ncommonnodes, &nparts,
        tpwgts, &ubvec, options,
        &edgecut, part, &comm);

    if (ret != METIS_OK) { /* handle error */ }

    MPI_Finalize();
    return 0;
}

Build docs developers (and LLMs) love