NODEND#
purpose#
The nodend
package finds a symmetric row and column permutation,
\(P A P^T\), of a symmetric, sparse matrix \(A\) with the aim of limiting
the fill-in during subsequent Cholesky-like factorization. The package
is actually a wrapper to the METIS_NodeND procedure from versions
4.0, 5.1 and 5.2 of the METIS package from the Karypis Lab.
See Section 4 of $GALAHAD/doc/nodend.pdf for additional details.
method#
Variants of node-based nested-dissection ordering are employed.
reference#
The basic methods are those described in
G. Karypis. METIS, A software package for partitioning unstructured graphs, partitioning meshes, and computing fill-reducing orderings of sparse matrices, Version 5, Department of Computer Science & Engineering, University of Minnesota Minneapolis, MN 55455, USA (2013), see
matrix storage#
The sparsity pattern of the symmetric \(n\) by \(n\) sparse matrix \(A\) may be presented and stored in a couple of formats. But crucially symmetry is exploited by only indices from the lower triangular part (i.e, those entries that lie on or below the leading diagonal).
Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 \leq l \leq ne-1\), of \(H\), its row index i and column index j, \(0 \leq j \leq i \leq n-1\), are stored as the \(l\)-th components of the integer arrays H_row and H_col, respectively, while the number of nonzeros is recorded as H_ne = \(ne\). Note that only the entries in the lower triangle should be stored. The string H_type = ‘coordinate’ should be specified.
Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(H\) the i-th component of the integer array H_ptr holds the position of the first entry in this row, while H_ptr(n) holds the total number of entries. The column indices j, \(0 \leq j \leq i\), are stored in components l = H_ptr(i), …, H_ptr(i+1)-1 of the integer array H_col. Note that as before only the entries in the lower triangle should be stored. For sparse matrices, this scheme almost always requires less storage than its predecessor. The string H_type = ‘sparse_by_rows’ should be specified.
functions#
- nodend.initialize()#
Set default option values and initialize private data.
Returns:
- optionsdict
- dictionary containing default control options:
- versionstr
specify the version of METIS to be used. Possible values are 4.0, 5.1 and 5.2.
- errorint
error and warning diagnostics occur on stream error.
- outint
general output occurs on stream out.
- print_levelint
the level of output required is specified by print_level.
- no_metis_4_use_5_insteadbool
if `` no_metis_4_use_5_instead`` is True, and METIS 4.0 is not availble, use Metis 5.2 instead.
- prefixstr
all output lines will be prefixed by the string contained in quotes within
prefix
, e.g. ‘word’ (note the qutoes) will result in the prefix word.- metis4_ptypeint
- the partitioning method employed. 0 = multilevel recursive
bisectioning: 1 = multilevel k-way partitioning
- metis4_ctypeint
the matching scheme to be used during coarsening: 1 = random matching, 2 = heavy-edge matching, 3 = sorted heavy-edge matching, and 4 = k-way sorted heavy-edge matching.
- metis4_itypeint
the algorithm used during initial partitioning: 1 = edge-based region growing and 2 = node-based region growing.
- metis4_rtypeint
the algorithm used for refinement: 1 = two-sided node Fiduccia-Mattheyses (FM) refinement, and 2 = one-sided node FM refinement.
- metis4_dbglvlint
the amount of progress/debugging information printed: 0 = nothing, 1 = timings, and \(>\) 1 increasingly more.
- metis4_oflagsint
select whether or not to compress the graph, and to order connected components separately: 0 = do neither, 1 = try to compress the graph, 2 = order each connected component separately, and 3 = do both.
- metis4_pfactorint
the minimum degree of the vertices that will be ordered last. More specifically, any vertices with a degree greater than 0.1 metis4_pfactor times the average degree are removed from the graph, an ordering of the rest of the vertices is computed, and an overall ordering is computed by ordering the removed vertices at the end of the overall ordering. Any value smaller than 1 means that no vertices will be ordered last.
- metis4_nsepsint
the number of different separators that the algorithm will compute at each level of nested dissection.
- metis5_ptypeint
the partitioning method. The value 0 gives multilevel recursive bisectioning, while 1 corresponds to multilevel \(k\)-way partitioning.
- metis5_objtypeint
the type of the objective. Currently the only and default value metis5_objtype = 2, specifies node-based nested dissection, and any invalid value will be replaced by this default.
- metis5_ctypeint
- the matching scheme to be used during coarsening:
0 = random matching, and 1 = sorted heavy-edge matching.
- metis5_iptypeint
the algorithm used during initial partitioning: 2 = derive separators from edge cuts, and 3 = grow bisections using a greedy node-based strategy.
- metis5_rtypeint
- the algorithm used for refinement: 2 = Two-sided node FM refinement,
and 3 = One-sided node FM refinement.
- metis5_dbglvlint
the amount of progress/debugging information printed: 0 = nothing, 1 = diagnostics, 2 = plus timings, and \(>\) 2 plus more.
- metis5_nipartsint
the number of initial partitions used by MeTiS 5.2.
- metis5_niterint
the number of iterations used by the refinement algorithm.
- metis5_ncutsint
the number of different partitionings that it will compute: -1 = not used.
- metis5_seedint
the seed for the random number generator.
- metis5_ondiskint
whether on-disk storage is used (0 = no, 1 = yes) by MeTiS 5.2.
- metis5_minconnint
specify that the partitioning routines should try to minimize the maximum degree of the subdomain graph: 0 = no, 1 = yes, and -1 = not used.
- metis5_contigint
specify that the partitioning routines should try to produce partitions that are contiguous: 0 = no, 1 = yes, and -1 = not used.
- metis5_compressint
specify that the graph should be compressed by combining together vertices that have identical adjacency lists: 0 = no, and 1 = yes.
- metis5_ccorderint
specify if the connected components of the graph should first be identified and ordered separately: 0 = no, and 1 = yes.
- metis5_pfactorint
the minimum degree of the vertices that will be ordered last. More specifically, any vertices with a degree greater than 0.1 metis4_pfactor times the average degree are removed from the graph, an ordering of the rest of the vertices is computed, and an overall ordering is computed by ordering the removed vertices at the end of the overall ordering.
- metis5_nsepsint
the number of different separators that the algorithm will compute at each level of nested dissection.
- metis5_ufactorint
the maximum allowed load imbalance (1 + metis5_ufactor)/1000 among the partitions.
- metis5_dropedgesint
will edges be dropped (0 = no, 1 = yes) by MeTiS 5.2.
- metis5_no2hopint
specify that the coarsening will not perform any 2–hop matchings when the standard matching approach fails to sufficiently coarsen the graph: 0 = no, and 1 = yes.
- metis5_twohopint
reserved for future use but ignored at present.
- metis5_fastint
reserved for future use but ignored at present.
- nodend.order(n, A_type, A_ne, A_row, A_col, A_ptr, options=None)#
Find a row/colum permutation for sparse Cholesky-like factorization.
Parameters:
- nint
holds the number of variables.
- A_typestring
specifies the symmetric storage scheme used for \(A\). It should be one of ‘coordinate’ or ‘sparse_by_rows’; lower or upper case variants are allowed. If A_type is not one of the supported values, the identity permutation will be returned.
- A_neint
holds the number of entries in the lower triangular part of \(A\) in the sparse co-ordinate storage scheme. It need not be set for the sparse-by-row scheme.
- A_rowndarray(A_ne)
holds the row indices of the lower triangular part of \(A\) in the sparse co-ordinate storage scheme. It need not be set for the sparse-by-row scheme, and in this case can be None.
- A_colndarray(A_ne)
holds the column indices of the lower triangular part of \(A\) in either of the supported storage schemes.
- A_ptrndarray(n+1)
holds the starting position of each row of the lower triangular part of \(A\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the coordinate scheme is used, and in this case can be None.
- optionsdict, optional
dictionary of control options (see
nodend.initialize
).Returns:
- permndarray(n)
holds the permutation array, so that the perm[\(i\)]-th rows and columns in the permuted matrix \(P A P^T\) correspond to those labelled i in \(A\), 0 \(\leq\) i \(\leq\) n-1.
- [optional] nodend.information()
Provide optional output information.
Returns:
- informdict
- dictionary containing output information:
- statusint
the return status. Possible values are:
0
The call was successful.
-1
An allocation error occurred. A message indicating the offending array is written on unit options[‘error’], and the returned allocation status and a string containing the name of the offending array are held in inform[‘alloc_status’] and inform[‘bad_alloc’] respectively.
-2
A deallocation error occurred. A message indicating the offending array is written on unit options[‘error’] and the returned allocation status and a string containing the name of the offending array are held in inform[‘alloc_status’] and inform[‘bad_alloc’] respectively.
-3
One of the restrictions n \(> 0\), A.n \(> 0\) or A.ne \(< 0\), for co-ordinate entry, or requirements that A.type contain its relevant string ‘coordinate’ or ‘sparse_by_rows’, and options[‘version’] in one of ‘4.0’, ‘5.1’ or ‘5.2’ has been violated.
-26
The requested version of METIS is not available.
-57
METIS has insufficient memory to continue.
-71
An internal METIS error occurred.
- alloc_statusint
the status of the last attempted allocation/deallocation.
- bad_allocstr
the name of the array for which an allocation/deallocation error occurred.
- versionstr
specifies the version of METIS that was actually used.
example code#
from galahad import nodend
import numpy as np
np.set_printoptions(precision=4,suppress=True,floatmode='fixed')
print("\n** python test: nodend")
# set parameters
n = 3
A_type = 'coordinate'
A_ne = 4
A_row = np.array([0,1,2,2])
A_col = np.array([0,1,1,2])
A_ptr = None
# allocate internal data and set default options
options = nodend.initialize()
# set some non-default options
options['print_level'] = 0
#print("options:", options)
# find the ordering
perm = nodend.order(n, A_type, A_ne, A_row, A_col, A_ptr, options)
print(" perm",perm)
# get information
inform = nodend.information()
print('** nodend exit status:', inform['status'])
This example code is available in $GALAHAD/src/nodend/Python/test_nodend.py .