petsc-3.12.5 2020-03-29
MatGetOrdering
Gets a reordering for a matrix to reduce fill or to improve numerical stability of LU factorization.
Synopsis
#include "petscmat.h"
PetscErrorCode MatGetOrdering(Mat mat,MatOrderingType type,IS *rperm,IS *cperm)
Collective on Mat
Input Parameters
| mat | - the matrix
|
| type | - type of reordering, one of the following:
|
MATORDERINGNATURAL - Natural
MATORDERINGND - Nested Dissection
MATORDERING1WD - One-way Dissection
MATORDERINGRCM - Reverse Cuthill-McKee
MATORDERINGQMD - Quotient Minimum Degree
Output Parameters
| rperm | - row permutation indices
|
| cperm | - column permutation indices
|
Options Database Key
-mat_view_ordering draw -plots matrix nonzero structure in new ordering
Notes
This DOES NOT actually reorder the matrix; it merely returns two index sets
that define a reordering. This is usually not used directly, rather use the
options PCFactorSetMatOrderingType()
The user can define additional orderings; see MatOrderingRegister().
These are generally only implemented for sequential sparse matrices.
The external packages that PETSc can use for direct factorization such as SuperLU do not accept orderings provided by
this call.
fill, reordering, natural, Nested Dissection,
One-way Dissection, Cholesky, Reverse Cuthill-McKee,
Quotient Minimum Degree
See Also
MatOrderingRegister(), PCFactorSetMatOrderingType()
Level
intermediate
Location
src/mat/order/sorder.c
Examples
src/mat/examples/tutorials/ex1.c.html
src/ksp/ksp/examples/tutorials/ex18.c.html
Implementations
MatGetOrdering_Flow_SeqAIJ in src/mat/impls/aij/seq/aijfact.c
Index of all MatOrderings routines
Table of Contents for all manual pages
Index of all manual pages