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_Mesh2Dual constructs the dual graph of a distributed finite element mesh in parallel. In the dual graph, each element becomes a vertex, and two vertices are connected by an edge when the corresponding elements share at least ncommonnodes mesh nodes — for example, 2 for edge-adjacent 2D elements or 3 for face-adjacent 3D elements. The resulting dual graph can then be passed to ParMETIS_V3_PartKway to partition the mesh.
The output arrays xadj and adjncy are allocated by ParMETIS. The caller is responsible for freeing them with free() after use.
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);

Parameters

elmdist
idx_t*
Distribution of elements across processes. elmdist[i] is the global index of the first element on process i, and elmdist[i+1] - elmdist[i] is the number of local elements. 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. Analogous to xadj in CSR graph format. 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].
numflag
idx_t*
Indexing convention. 0 for C-style 0-based indexing; 1 for Fortran-style 1-based indexing.
ncommonnodes
idx_t*
Minimum number of shared mesh nodes required for two elements to be considered adjacent in the dual graph. Use 2 for edge adjacency in 2D meshes and 3 for face adjacency in 3D meshes.
xadj
idx_t**
Output. Pointer to the allocated adjacency row pointer array of the dual graph, in CSR format. ParMETIS allocates this array; the caller must free it with free(*xadj).
adjncy
idx_t**
Output. Pointer to the allocated adjacency column index array of the dual graph, in CSR format. ParMETIS allocates this array; the caller must free it with free(*adjncy).
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 <stdlib.h>
#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;

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

    idx_t numflag       = 0;
    idx_t ncommonnodes  = 3; /* face adjacency for 3D */

    idx_t *dual_xadj   = NULL;
    idx_t *dual_adjncy = NULL;

    int ret = ParMETIS_V3_Mesh2Dual(
        elmdist, eptr, eind,
        &numflag, &ncommonnodes,
        &dual_xadj, &dual_adjncy,
        &comm);

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

    /* Use dual_xadj and dual_adjncy for partitioning ... */

    free(dual_xadj);
    free(dual_adjncy);

    MPI_Finalize();
    return 0;
}

Build docs developers (and LLMs) love