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_NodeND computes a fill-reducing permutation for a sparse symmetric matrix using the multilevel nested dissection algorithm. The resulting permutation, when applied to reorder the rows and columns of the matrix before factorization, significantly reduces the number of nonzeros introduced during Cholesky or LU decomposition (fill). It is the standard preprocessing step for direct sparse solvers such as CHOLMOD, PARDISO, SuperLU, and similar libraries.
The vwgt parameter is NULL in most practical use cases. Vertex weights are only needed when the graph has vertices representing unequal amounts of work in the downstream computation (e.g., supernodes of different sizes). When in doubt, pass NULL.

METIS_NodeND

Computes a fill-reducing ordering of the vertices of a graph. The graph should represent the nonzero structure of the lower triangle of the symmetric matrix to be factored (diagonal entries are excluded).
int METIS_NodeND(
    idx_t  *nvtxs,
    idx_t  *xadj,
    idx_t  *adjncy,
    idx_t  *vwgt,
    idx_t  *options,
    idx_t  *perm,
    idx_t  *iperm
);
nvtxs
idx_t *
required
Number of vertices in the graph (equivalently, the order of the sparse matrix).
xadj
idx_t *
required
CSR row-pointer array of size nvtxs + 1. xadj[i] is the index into adjncy where vertex i’s adjacency list starts.
adjncy
idx_t *
required
CSR column-index array. Stores all adjacency lists concatenated. The adjacency list of each vertex should not contain the vertex itself (no self-loops).
vwgt
idx_t *
Vertex weights array of size nvtxs. Controls how separator size is measured during nested dissection. Pass NULL in most cases; METIS will treat all vertices as having unit weight.
options
idx_t *
required
Options array of size METIS_NOPTIONS initialized with METIS_SetDefaultOptions. Key ordering-specific options include METIS_OPTION_NSEPS, METIS_OPTION_COMPRESS, METIS_OPTION_CCORDER, and METIS_OPTION_PFACTOR. See Options Reference.
perm
idx_t *
required
Output: permutation array of size nvtxs. If A is the original matrix and A' is the permuted matrix, then A'[i][j] = A[perm[i]][perm[j]]. Equivalently, perm[i] gives the original row/column index that maps to position i in the reordered matrix. Must be pre-allocated by the caller.
iperm
idx_t *
required
Output: inverse permutation array of size nvtxs. iperm[i] gives the new position of original vertex i in the reordered matrix. The relationship between the two arrays is: perm[iperm[i]] == i and iperm[perm[i]] == i. Must be pre-allocated by the caller.

Understanding perm and iperm

The two output arrays encode the same reordering from opposite directions:
perm[new_index]  = old_index   (new → old)
iperm[old_index] = new_index   (old → new)
Sparse solvers typically need one or the other depending on their convention. For example:
  • To reorder a matrix A to A' = P*A*P^T, use perm to index into the original row/column arrays.
  • To map a right-hand side vector b to b' = P*b, apply b'[iperm[i]] = b[i] for all i.

Using the Ordering with a Sparse Solver

1

Build the graph from your matrix

Extract the nonzero pattern of your symmetric matrix as a CSR graph. Diagonal entries and the upper triangle are excluded; only store the off-diagonal entries of the lower triangle (or use the full symmetric pattern — METIS handles both).
/* Matrix sparsity pattern -> CSR graph */
idx_t nvtxs = n;               /* matrix order */
idx_t *xadj   = /* row pointers, size n+1 */;
idx_t *adjncy = /* column indices, size nnz */;
2

Compute the fill-reducing ordering

idx_t perm[n], iperm[n];
idx_t options[METIS_NOPTIONS];
METIS_SetDefaultOptions(options);

int ret = METIS_NodeND(&nvtxs, xadj, adjncy, NULL, options, perm, iperm);
if (ret != METIS_OK) { /* handle error */ }
3

Reorder rows and columns

Apply the permutation to build the reordered matrix A' = P*A*P^T:
/* For each nonzero A[row][col], place it at A'[iperm[row]][iperm[col]] */
for (int i = 0; i < n; i++) {
    int new_row = iperm[i];
    for (int j = xadj[i]; j < xadj[i+1]; j++) {
        int new_col = iperm[adjncy[j]];
        /* insert A_values[j] into new_row, new_col of A' */
    }
}
Many sparse solver libraries accept the permutation arrays directly and perform the reordering internally.
4

Factor and solve

Pass the reordered matrix to your sparse solver. After solving the permuted system A'*x' = P*b, recover the original solution with x[perm[i]] = x'[i].

Complete C Example

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

int main(void) {
    /*
     * Symmetric matrix sparsity pattern (off-diagonal, lower triangle):
     *
     *   [x . x . .]      row 0: col 2
     *   [. x x . .]      row 1: col 2
     *   [x x x x .]      row 2: cols 0, 1, 3
     *   [. . x x x]      row 3: cols 2, 4
     *   [. . . x x]      row 4: col 3
     *
     * As undirected graph (both directions):
     */
    idx_t nvtxs = 5;
    idx_t xadj[]   = {0, 1, 2, 5, 8, 10};
    idx_t adjncy[] = {
        2,        /* vertex 0: neighbor 2 */
        2,        /* vertex 1: neighbor 2 */
        0, 1, 3,  /* vertex 2: neighbors 0, 1, 3 */
        2, 4,     /* vertex 3: neighbors 2, 4 */
        3         /* vertex 4: neighbor 3 */
        /* (symmetric: 5 undirected edges stored as 10 directed entries) */
    };

    /* Adjust adjncy to include both directions for METIS */
    idx_t xadj_sym[]   = {0, 1, 2, 5, 8, 10};
    idx_t adjncy_sym[] = {2, 2, 0,1,3, 2,4, 3};

    idx_t perm[5], iperm[5];
    idx_t options[METIS_NOPTIONS];
    METIS_SetDefaultOptions(options);

    int ret = METIS_NodeND(&nvtxs, xadj_sym, adjncy_sym,
                           NULL,   /* unit vertex weights */
                           options,
                           perm, iperm);

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

    printf("Permutation (new -> old):\n");
    for (int i = 0; i < nvtxs; i++)
        printf("  perm[%d] = %d\n", i, (int)perm[i]);

    printf("Inverse permutation (old -> new):\n");
    for (int i = 0; i < nvtxs; i++)
        printf("  iperm[%d] = %d\n", i, (int)iperm[i]);

    return 0;
}

METIS_NodeNDP

A variant of METIS_NodeND designed for use by ParMETIS. It computes fill-reducing orderings for parallel direct solvers by performing a nested dissection that is aware of the desired number of processors. In addition to the permutation and inverse permutation, it returns separator sizes at each level of the dissection tree.
int METIS_NodeNDP(
    idx_t   nvtxs,
    idx_t  *xadj,
    idx_t  *adjncy,
    idx_t  *vwgt,
    idx_t   npes,
    idx_t  *options,
    idx_t  *perm,
    idx_t  *iperm,
    idx_t  *sizes
);
nvtxs
idx_t
required
Number of vertices. Passed by value (not a pointer), unlike METIS_NodeND.
xadj
idx_t *
required
CSR row-pointer array of size nvtxs + 1.
adjncy
idx_t *
required
CSR column-index array.
vwgt
idx_t *
Vertex weights array of size nvtxs. Pass NULL for unit weights.
npes
idx_t
required
Number of parallel processors. The dissection tree will have 2*npes - 1 nodes. Passed by value.
options
idx_t *
required
Options array initialized with METIS_SetDefaultOptions.
perm
idx_t *
required
Output: permutation array of size nvtxs. Pre-allocate before calling.
iperm
idx_t *
required
Output: inverse permutation array of size nvtxs. Pre-allocate before calling.
sizes
idx_t *
required
Output: array of size 2*npes that stores the sizes of the vertex ranges for each subgraph and separator in the dissection tree. Used by ParMETIS to determine communication patterns. Pre-allocate before calling.
METIS_NodeNDP is designed as an internal utility for ParMETIS and direct parallel solver interfaces. For standalone fill-reducing ordering, use METIS_NodeND.

Build docs developers (and LLMs) love