xref: /petsc/include/petscmatcoarsen.h (revision bae903cb6e5fe998dd743ea97ba6d0e4e52494ec)
126bd1501SBarry Smith #ifndef PETSCMATCOARSEN_H
226bd1501SBarry Smith #define PETSCMATCOARSEN_H
3ac09b921SBarry Smith 
4484f0a72SBarry Smith #include <petscmat.h>
5484f0a72SBarry Smith 
6ac09b921SBarry Smith /* SUBMANSEC = Mat */
7ac09b921SBarry Smith 
8484f0a72SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList;
9484f0a72SBarry Smith 
10484f0a72SBarry Smith /*S
11484f0a72SBarry Smith      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
12484f0a72SBarry Smith 
13484f0a72SBarry Smith    Level: advanced
14484f0a72SBarry Smith 
1595452b02SPatrick Sanan   Notes:
1695452b02SPatrick Sanan     This is used by the PCGAMG to generate coarser representations of an algebraic problem
17484f0a72SBarry Smith 
18db781477SPatrick Sanan .seealso: `MatCoarsenCreate()`, `MatCoarsenType`
19484f0a72SBarry Smith S*/
20484f0a72SBarry Smith typedef struct _p_MatCoarsen* MatCoarsen;
21484f0a72SBarry Smith 
22484f0a72SBarry Smith /*J
23484f0a72SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen algorithm
24484f0a72SBarry Smith 
25484f0a72SBarry Smith    Level: beginner
26484f0a72SBarry Smith 
27db781477SPatrick Sanan .seealso: `MatCoarsenCreate()`, `MatCoarsen`
28484f0a72SBarry Smith J*/
29484f0a72SBarry Smith typedef const char* MatCoarsenType;
30484f0a72SBarry Smith #define MATCOARSENMIS  "mis"
31484f0a72SBarry Smith #define MATCOARSENHEM  "hem"
32*bae903cbSmarkadams4 #define MATCOARSENMISK "misk"
33484f0a72SBarry Smith 
34484f0a72SBarry Smith /* linked list for aggregates */
35484f0a72SBarry Smith typedef struct _PetscCDIntNd{
36484f0a72SBarry Smith   struct _PetscCDIntNd *next;
37484f0a72SBarry Smith   PetscInt             gid;
38484f0a72SBarry Smith }PetscCDIntNd;
39484f0a72SBarry Smith 
40484f0a72SBarry Smith /* only used by node pool */
41484f0a72SBarry Smith typedef struct _PetscCDArrNd{
42484f0a72SBarry Smith   struct _PetscCDArrNd *next;
43484f0a72SBarry Smith   struct _PetscCDIntNd *array;
44484f0a72SBarry Smith }PetscCDArrNd;
45484f0a72SBarry Smith 
46*bae903cbSmarkadams4 /* 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 */
47484f0a72SBarry Smith typedef struct _PetscCoarsenData{
48484f0a72SBarry Smith   PetscCDArrNd pool_list;  /* node pool */
49484f0a72SBarry Smith   PetscCDIntNd *new_node;
50484f0a72SBarry Smith   PetscInt     new_left;
51*bae903cbSmarkadams4   PetscInt     chk_sz; /* chunck size */
52484f0a72SBarry Smith   PetscCDIntNd *extra_nodes;
53484f0a72SBarry Smith   PetscCDIntNd **array;  /* Array of lists */
54*bae903cbSmarkadams4   PetscInt     size; /* size of 'array' */
55484f0a72SBarry Smith   Mat          mat;  /* cache a Mat for communication data */
56484f0a72SBarry Smith }PetscCoarsenData;
57484f0a72SBarry Smith 
58484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
59484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
60484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
61484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
62484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
63484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetData(MatCoarsen, PetscCoarsenData **);
64484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
65484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
66484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
67484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
68484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
69484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
70fe2efc57SMark PETSC_EXTERN PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen,PetscObject,const char[]);
71484f0a72SBarry Smith 
72484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt,PetscCoarsenData**);
73484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData*);
74484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd*,PetscInt);
75484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd*,PetscInt*);
76484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData*,PetscInt,PetscInt);
77484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDAppendRemove(PetscCoarsenData*,PetscInt,PetscInt);
78484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData*,PetscInt,PetscCDIntNd*);
79484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData*,PetscInt,PetscCDIntNd*);
80484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDSizeAt(const PetscCoarsenData*,PetscInt,PetscInt*);
81484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDEmptyAt(const PetscCoarsenData*,PetscInt,PetscBool*);
82484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDSetChuckSize(PetscCoarsenData*,PetscInt);
83484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData*,MPI_Comm);
84484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDGetMIS(PetscCoarsenData*,IS*);
85*bae903cbSmarkadams4 PETSC_EXTERN PetscErrorCode PetscCDGetMat(PetscCoarsenData*,Mat*);
86484f0a72SBarry Smith PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData*,Mat);
87*bae903cbSmarkadams4 PETSC_EXTERN PetscErrorCode PetscCDRemoveAll(PetscCoarsenData*, PetscInt);
88484f0a72SBarry Smith 
89539c167fSBarry Smith PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData*,PetscInt,PetscCDIntNd**);
90539c167fSBarry Smith PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData*,PetscInt,PetscCDIntNd**);
910a3c815dSMark Adams PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData*,const PetscInt,Mat,PetscInt*,IS**);
92484f0a72SBarry Smith 
93484f0a72SBarry Smith #endif
94