Actual source code: spartition.c
2: #include <petscmat.h>
3: #include <petsc/private/matimpl.h>
5: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Current(MatPartitioning);
6: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Average(MatPartitioning part);
7: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Square(MatPartitioning);
8: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Parmetis(MatPartitioning);
9: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning);
10: #if defined(PETSC_HAVE_CHACO)
11: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Chaco(MatPartitioning);
12: #endif
13: #if defined(PETSC_HAVE_PARTY)
14: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning);
15: #endif
16: #if defined(PETSC_HAVE_PTSCOTCH)
17: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_PTScotch(MatPartitioning);
18: #endif
20: /*@C
21: MatPartitioningRegisterAll - Registers all of the matrix partitioning routines in PETSc.
23: Not Collective
25: Level: developer
27: Adding new methods:
28: To add a new method to the registry. Copy this routine and
29: modify it to incorporate a call to `MatPartitioningRegister()` for
30: the new method, after the current list.
32: Restricting the choices: To prevent all of the methods from being
33: registered and thus save memory, copy this routine and modify it to
34: register a zero, instead of the function name, for those methods you
35: do not wish to register. Make sure that the replacement routine is
36: linked before libpetscmat.a.
38: .seealso: `MatPartitioning`, `MatPartitioningType`, `MatPartitioningRegister()`, `MatPartitioningRegisterDestroy()`
39: @*/
40: PetscErrorCode MatPartitioningRegisterAll(void)
41: {
42: PetscFunctionBegin;
43: if (MatPartitioningRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
44: MatPartitioningRegisterAllCalled = PETSC_TRUE;
46: PetscCall(MatPartitioningRegister(MATPARTITIONINGCURRENT, MatPartitioningCreate_Current));
47: PetscCall(MatPartitioningRegister(MATPARTITIONINGAVERAGE, MatPartitioningCreate_Average));
48: PetscCall(MatPartitioningRegister(MATPARTITIONINGSQUARE, MatPartitioningCreate_Square));
49: PetscCall(MatPartitioningRegister(MATPARTITIONINGHIERARCH, MatPartitioningCreate_Hierarchical));
50: #if defined(PETSC_HAVE_PARMETIS)
51: PetscCall(MatPartitioningRegister(MATPARTITIONINGPARMETIS, MatPartitioningCreate_Parmetis));
52: #endif
53: #if defined(PETSC_HAVE_CHACO)
54: PetscCall(MatPartitioningRegister(MATPARTITIONINGCHACO, MatPartitioningCreate_Chaco));
55: #endif
56: #if defined(PETSC_HAVE_PARTY)
57: PetscCall(MatPartitioningRegister(MATPARTITIONINGPARTY, MatPartitioningCreate_Party));
58: #endif
59: #if defined(PETSC_HAVE_PTSCOTCH)
60: PetscCall(MatPartitioningRegister(MATPARTITIONINGPTSCOTCH, MatPartitioningCreate_PTScotch));
61: #endif
62: PetscFunctionReturn(PETSC_SUCCESS);
63: }