xref: /petsc/include/petscmatcoarsen.h (revision df4cd43f92eaa320656440c40edb1046daee8f75)
1 #ifndef PETSCMATCOARSEN_H
2 #define PETSCMATCOARSEN_H
3 
4 #include <petscmat.h>
5 
6 /* SUBMANSEC = Mat */
7 
8 PETSC_EXTERN PetscFunctionList MatCoarsenList;
9 
10 /*S
11      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
12 
13    Level: advanced
14 
15   Note:
16     This is used by the `PCGAMG` to generate coarser representations of an algebraic problem
17 
18 .seealso: [](chapter_matrices), [](sec_graph), `Mat`, `MatCoarsenCreate()`, `MatCoarsenType`, `MatColoringType`, `MatPartitioningType`, `MatOrderingType`
19           `MatColoring`, `MatPartitioning`
20 S*/
21 typedef struct _p_MatCoarsen *MatCoarsen;
22 
23 /*J
24     MatCoarsenType - String with the name of a PETSc matrix coarsening algorithm
25 
26    Level: beginner
27 
28 .seealso: [](chapter_matrices), [](sec_graph), `Mat`, `MatCoarsenCreate()`, `MatCoarsen`, `MatColoringType`, `MatPartitioningType`, `MatOrderingType`
29 J*/
30 typedef const char *MatCoarsenType;
31 #define MATCOARSENMIS  "mis"
32 #define MATCOARSENHEM  "hem"
33 #define MATCOARSENMISK "misk"
34 
35 /* linked list for aggregates */
36 typedef struct _PetscCDIntNd {
37   struct _PetscCDIntNd *next;
38   PetscInt              gid;
39 } PetscCDIntNd;
40 
41 /* only used by node pool */
42 typedef struct _PetscCDArrNd {
43   struct _PetscCDArrNd *next;
44   struct _PetscCDIntNd *array;
45 } PetscCDArrNd;
46 
47 /* linked list data structure that encodes aggragates and C-F points with array[idx] == NULL for F point and array of indices in an aggrate or C point (first index is always global index my0 + idx */
48 typedef struct _PetscCoarsenData {
49   PetscCDArrNd   pool_list; /* node pool */
50   PetscCDIntNd  *new_node;
51   PetscInt       new_left;
52   PetscInt       chk_sz; /* chunk size */
53   PetscCDIntNd  *extra_nodes;
54   PetscCDIntNd **array; /* Array of lists */
55   PetscInt       size;  /* size of 'array' */
56   Mat            mat;   /* cache a Mat for communication data */
57 } PetscCoarsenData;
58 
59 PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm, MatCoarsen *);
60 PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen, MatCoarsenType);
61 PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen, Mat);
62 PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen, const IS);
63 PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen, PetscBool);
64 PETSC_EXTERN PetscErrorCode MatCoarsenGetData(MatCoarsen, PetscCoarsenData **);
65 PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
66 PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen *);
67 PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[], PetscErrorCode (*)(MatCoarsen));
68 PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen, PetscViewer);
69 PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
70 PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen, MatCoarsenType *);
71 PETSC_EXTERN PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen, PetscObject, const char[]);
72 
73 PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt, PetscCoarsenData **);
74 PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData *);
75 PETSC_EXTERN PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *, PetscInt);
76 PETSC_EXTERN PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *, PetscInt *);
77 PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData *, PetscInt, PetscInt);
78 PETSC_EXTERN PetscErrorCode PetscCDAppendRemove(PetscCoarsenData *, PetscInt, PetscInt);
79 PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
80 PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
81 PETSC_EXTERN PetscErrorCode PetscCDSizeAt(const PetscCoarsenData *, PetscInt, PetscInt *);
82 PETSC_EXTERN PetscErrorCode PetscCDEmptyAt(const PetscCoarsenData *, PetscInt, PetscBool *);
83 PETSC_EXTERN PetscErrorCode PetscCDSetChuckSize(PetscCoarsenData *, PetscInt);
84 PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData *, MPI_Comm);
85 PETSC_EXTERN PetscErrorCode PetscCDGetMIS(PetscCoarsenData *, IS *);
86 PETSC_EXTERN PetscErrorCode PetscCDGetMat(PetscCoarsenData *, Mat *);
87 PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData *, Mat);
88 PETSC_EXTERN PetscErrorCode PetscCDRemoveAll(PetscCoarsenData *, PetscInt);
89 
90 PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
91 PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
92 PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, Mat, PetscInt *, IS **);
93 
94 PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
95 PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);
96 #endif
97