Actual source code: spartition.c
1: #include <petscmat.h>
2: #include <petsc/private/matimpl.h>
4: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Current(MatPartitioning);
5: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Average(MatPartitioning part);
6: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Square(MatPartitioning);
7: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Parmetis(MatPartitioning);
8: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning);
9: #if defined(PETSC_HAVE_CHACO)
10: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Chaco(MatPartitioning);
11: #endif
12: #if defined(PETSC_HAVE_PARTY)
13: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning);
14: #endif
15: #if defined(PETSC_HAVE_PTSCOTCH)
16: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_PTScotch(MatPartitioning);
17: #endif
19: /*@C
20: MatPartitioningRegisterAll - Registers all of the matrix partitioning routines in PETSc.
22: Not Collective
24: Level: developer
26: .seealso: `MatPartitioning`, `MatPartitioningType`, `MatPartitioningRegister()`, `MatPartitioningRegisterDestroy()`
27: @*/
28: PetscErrorCode MatPartitioningRegisterAll(void)
29: {
30: PetscFunctionBegin;
31: if (MatPartitioningRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
32: MatPartitioningRegisterAllCalled = PETSC_TRUE;
34: PetscCall(MatPartitioningRegister(MATPARTITIONINGCURRENT, MatPartitioningCreate_Current));
35: PetscCall(MatPartitioningRegister(MATPARTITIONINGAVERAGE, MatPartitioningCreate_Average));
36: PetscCall(MatPartitioningRegister(MATPARTITIONINGSQUARE, MatPartitioningCreate_Square));
37: PetscCall(MatPartitioningRegister(MATPARTITIONINGHIERARCH, MatPartitioningCreate_Hierarchical));
38: #if defined(PETSC_HAVE_PARMETIS)
39: PetscCall(MatPartitioningRegister(MATPARTITIONINGPARMETIS, MatPartitioningCreate_Parmetis));
40: #endif
41: #if defined(PETSC_HAVE_CHACO)
42: PetscCall(MatPartitioningRegister(MATPARTITIONINGCHACO, MatPartitioningCreate_Chaco));
43: #endif
44: #if defined(PETSC_HAVE_PARTY)
45: PetscCall(MatPartitioningRegister(MATPARTITIONINGPARTY, MatPartitioningCreate_Party));
46: #endif
47: #if defined(PETSC_HAVE_PTSCOTCH)
48: PetscCall(MatPartitioningRegister(MATPARTITIONINGPTSCOTCH, MatPartitioningCreate_PTScotch));
49: #endif
50: PetscFunctionReturn(PETSC_SUCCESS);
51: }