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