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.
The METIS C API provides a small set of high-level functions for graph and mesh partitioning and for computing fill-reducing orderings of sparse matrices. All public functions share a consistent interface: every parameter is passed by pointer, an idx_t options[METIS_NOPTIONS] array controls algorithm behavior, and every function returns one of four integer status codes.
Including and Linking
Add metis.h to your source file and link against libmetis:
# Static library (default build)
gcc myapp.c -o myapp -lmetis
# If the library is installed in a non-standard location
gcc myapp.c -o myapp -I/path/to/include -L/path/to/lib -lmetis
idx_t and real_t are configured at library build time via the IDXTYPEWIDTH and REALTYPEWIDTH preprocessor defines. idx_t is either int32_t (32-bit) or int64_t (64-bit). real_t is either float or double. Your application code must be compiled against the same metis.h that was used to build the library, or the type widths will not match and results will be undefined.
Return Codes
Every METIS function returns one of the following integer values defined in the rstatus_et enum:
| Constant | Value | Meaning |
|---|
METIS_OK | 1 | Function completed successfully |
METIS_ERROR_INPUT | -2 | Invalid input parameters or options |
METIS_ERROR_MEMORY | -3 | Insufficient memory |
METIS_ERROR | -4 | Other internal error |
Always check the return value before using output arrays.
The Options Array
All partitioning and ordering functions accept an idx_t options[METIS_NOPTIONS] array, where METIS_NOPTIONS equals 40. The recommended pattern is to call METIS_SetDefaultOptions first, which fills the array with -1 (meaning “use default for each entry”), then override only the options you need.
idx_t options[METIS_NOPTIONS];
METIS_SetDefaultOptions(options);
/* Override specific options */
options[METIS_OPTION_NCUTS] = 3; /* try 3 different cuts, keep the best */
options[METIS_OPTION_SEED] = 42; /* reproducible run */
options[METIS_OPTION_CONTIG] = 1; /* require contiguous partitions */
Passing NULL as the options argument is not supported; always pass a valid array initialized by METIS_SetDefaultOptions. See the Options Array Reference for the full list of option constants and their valid values.
Minimal Working Example
#include <stdio.h>
#include <stdlib.h>
#include <metis.h>
int main(void) {
/* A simple 6-vertex graph in CSR format */
idx_t nvtxs = 6, ncon = 1, nparts = 2;
idx_t xadj[] = {0, 2, 5, 7, 9, 12, 14};
idx_t adjncy[] = {1,2, 0,2,3, 0,1, 1,4, 3,5,0, 4,3};
idx_t edgecut;
idx_t part[6];
idx_t options[METIS_NOPTIONS];
int ret;
METIS_SetDefaultOptions(options);
ret = METIS_PartGraphKway(
&nvtxs, &ncon,
xadj, adjncy,
NULL, /* vwgt — unit weights */
NULL, /* vsize — not used */
NULL, /* adjwgt — unit edge weights */
&nparts,
NULL, /* tpwgts — equal target weights */
NULL, /* ubvec — default imbalance */
options,
&edgecut, part
);
if (ret != METIS_OK) {
fprintf(stderr, "METIS partitioning failed with code %d\n", ret);
return 1;
}
printf("Edge cut: %d\n", (int)edgecut);
for (idx_t i = 0; i < nvtxs; i++)
printf(" vertex %d -> partition %d\n", (int)i, (int)part[i]);
return 0;
}
All graph-partitioning functions represent graphs in Compressed Sparse Row (CSR) format using two arrays:
xadj — integer array of size nvtxs + 1. xadj[i] is the index into adjncy where vertex i’s adjacency list begins; xadj[i+1] is where it ends. So vertex i has degree xadj[i+1] - xadj[i].
adjncy — integer array of total length equal to the sum of all vertex degrees. Stores the neighbor lists of all vertices back-to-back.
For an undirected graph, every edge (u, v) must appear twice: once in u’s adjacency list and once in v’s. Edge weights, when provided via adjwgt, follow the same layout as adjncy.
Vertex 0: neighbors at adjncy[xadj[0] .. xadj[1]-1]
Vertex 1: neighbors at adjncy[xadj[1] .. xadj[2]-1]
...
Fortran Bindings
METIS provides Fortran-callable wrappers for every public API function. Four naming variants are generated per function to accommodate different Fortran compilers:
| Variant | Example |
|---|
| All uppercase | METIS_PARTGRAPHKWAY |
| All lowercase | metis_partgraphkway |
| Lowercase with single trailing underscore | metis_partgraphkway_ |
| Lowercase with double trailing underscore | metis_partgraphkway__ |
Fortran callers pass all arguments by reference (pointer), matching the C API’s pointer-based convention. Index numbering should be set to 1-based via METIS_OPTION_NUMBERING = 1.
Public Function Summary
| Function | Description | Reference |
|---|
METIS_PartGraphKway | Partition a graph into k parts using multilevel k-way partitioning | Graph Partitioning |
METIS_PartGraphRecursive | Partition a graph using multilevel recursive bisection | Graph Partitioning |
METIS_PartMeshDual | Partition a mesh by partitioning its dual graph | Mesh Partitioning |
METIS_PartMeshNodal | Partition a mesh by partitioning its nodal graph | Mesh Partitioning |
METIS_MeshToDual | Convert a mesh to its dual graph | Mesh Partitioning |
METIS_MeshToNodal | Convert a mesh to its nodal graph | Mesh Partitioning |
METIS_NodeND | Compute a fill-reducing ordering using nested dissection | Ordering |
METIS_SetDefaultOptions | Initialize an options array with default values | Options |
METIS_Free | Free memory allocated and returned by METIS | Data Types |
METIS_ComputeVertexSeparator | Compute a vertex separator of a graph (used by ParMETIS) | Graph Partitioning |
METIS_NodeNDP | Nested dissection ordering for ParMETIS | Ordering |