xref: /petsc/src/mat/graphops/order/sregis.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9)
1*8be712e4SBarry Smith #include <petsc/private/matimpl.h> /*I       "petscmat.h"   I*/
2*8be712e4SBarry Smith 
3*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_Natural(Mat, MatOrderingType, IS *, IS *);
4*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_ND(Mat, MatOrderingType, IS *, IS *);
5*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_1WD(Mat, MatOrderingType, IS *, IS *);
6*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_QMD(Mat, MatOrderingType, IS *, IS *);
7*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_RCM(Mat, MatOrderingType, IS *, IS *);
8*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_RowLength(Mat, MatOrderingType, IS *, IS *);
9*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_DSC(Mat, MatOrderingType, IS *, IS *);
10*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_WBM(Mat, MatOrderingType, IS *, IS *);
11*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_Spectral(Mat, MatOrderingType, IS *, IS *);
12*8be712e4SBarry Smith #if defined(PETSC_HAVE_SUITESPARSE)
13*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_AMD(Mat, MatOrderingType, IS *, IS *);
14*8be712e4SBarry Smith #endif
15*8be712e4SBarry Smith #if defined(PETSC_HAVE_METIS)
16*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_METISND(Mat, MatOrderingType, IS *, IS *);
17*8be712e4SBarry Smith #endif
18*8be712e4SBarry Smith 
19*8be712e4SBarry Smith /*@C
20*8be712e4SBarry Smith   MatOrderingRegisterAll - Registers all of the matrix
21*8be712e4SBarry Smith   reordering routines in PETSc.
22*8be712e4SBarry Smith 
23*8be712e4SBarry Smith   Not Collective
24*8be712e4SBarry Smith 
25*8be712e4SBarry Smith   Level: developer
26*8be712e4SBarry Smith 
27*8be712e4SBarry Smith   Notes:
28*8be712e4SBarry Smith   To add a new method to the registry. Copy this routine and
29*8be712e4SBarry Smith   modify it to incorporate a call to `MatReorderRegister()` for
30*8be712e4SBarry Smith   the new method, after the current list.
31*8be712e4SBarry Smith 
32*8be712e4SBarry Smith   To prevent all of the methods from being
33*8be712e4SBarry Smith   registered and thus save memory, copy this routine and comment out
34*8be712e4SBarry Smith   those orderigs you do not wish to include.  Make sure that the
35*8be712e4SBarry Smith   replacement routine is linked before libpetscmat.a.
36*8be712e4SBarry Smith 
37*8be712e4SBarry Smith .seealso: `MatOrderingType`, `MatOrderingRegister()`
38*8be712e4SBarry Smith @*/
MatOrderingRegisterAll(void)39*8be712e4SBarry Smith PetscErrorCode MatOrderingRegisterAll(void)
40*8be712e4SBarry Smith {
41*8be712e4SBarry Smith   PetscFunctionBegin;
42*8be712e4SBarry Smith   if (MatOrderingRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
43*8be712e4SBarry Smith   MatOrderingRegisterAllCalled = PETSC_TRUE;
44*8be712e4SBarry Smith 
45*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGNATURAL, MatGetOrdering_Natural));
46*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGND, MatGetOrdering_ND));
47*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERING1WD, MatGetOrdering_1WD));
48*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGRCM, MatGetOrdering_RCM));
49*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGQMD, MatGetOrdering_QMD));
50*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGROWLENGTH, MatGetOrdering_RowLength));
51*8be712e4SBarry Smith #if defined(PETSC_HAVE_SUPERLU_DIST)
52*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGWBM, MatGetOrdering_WBM));
53*8be712e4SBarry Smith #endif
54*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGSPECTRAL, MatGetOrdering_Spectral));
55*8be712e4SBarry Smith #if defined(PETSC_HAVE_SUITESPARSE)
56*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGAMD, MatGetOrdering_AMD));
57*8be712e4SBarry Smith #endif
58*8be712e4SBarry Smith #if defined(PETSC_HAVE_METIS)
59*8be712e4SBarry Smith   PetscCall(MatOrderingRegister(MATORDERINGMETISND, MatGetOrdering_METISND));
60*8be712e4SBarry Smith #endif
61*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
62*8be712e4SBarry Smith }
63