xref: /petsc/include/petscdmforest.h (revision 163d334ea1df412efbe58c2e8fce360a2afeaf87)
1 /*
2   DMForest, for parallel, hierarchically refined, distributed mesh problems.
3 */
4 #if !defined(__PETSCDMFOREST_H)
5 #define __PETSCDMFOREST_H
6 
7 #include <petscdm.h>
8 
9 /*J
10     DMForestTopology - String with the name of a PETSc DMForest base mesh topology
11 
12    Level: beginner
13 
14 .seealso: DMForestSetType(), DMForest
15 J*/
16 typedef const char* DMForestTopology;
17 
18 /* Just a name for the shape of the domain */
19 PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology);
20 PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *);
21 
22 /* this is the coarsest possible forest: can be any DM which we can
23  * convert to a DMForest (right now: plex) */
24 PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM);
25 PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *);
26 
27 /* this is the forest from which we adapt */
28 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityForest(DM, DM);
29 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityForest(DM, DM *);
30 
31 /* reserve some adaptivity types */
32 enum {DM_FOREST_KEEP = 0,
33       DM_FOREST_REFINE,
34       DM_FOREST_COARSEN,
35       DM_FOREST_RESERVED_ADAPTIVITY_COUNT};
36 
37 typedef PetscInt DMForestAdaptivityPurpose;
38 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityPurpose(DM, DMForestAdaptivityPurpose);
39 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityPurpose(DM, DMForestAdaptivityPurpose*);
40 
41 /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */
42 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt);
43 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *);
44 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt);
45 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *);
46 
47 PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt);
48 PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *);
49 
50 PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt);
51 PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *);
52 
53 PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt);
54 PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *);
55 
56 PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt);
57 PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *);
58 
59 PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *);
60 PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *);
61 
62 
63 /* flag each cell with an adaptivity count: should match the cell section */
64 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityLabel(DM, const char *);
65 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityLabel(DM, const char **);
66 
67 /*J
68     DMForestAdaptivityStrategy - String with the name of a PETSc DMForest adaptivity strategy
69 
70    Level: intermediate
71 
72 .seealso: DMForestSetType(), DMForest
73 J*/
74 typedef const char* DMForestAdaptivityStrategy;
75 #define DMFORESTADAPTALL "all"
76 #define DMFORESTADAPTANY "any"
77 
78 /* how to combine: -flags         from multiple processes,
79  *                 -coarsen flags from multiple children
80  */
81 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy);
82 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *);
83 
84 PETSC_EXTERN PetscErrorCode DMForestSetComputeAdaptivitySF(DM, PetscBool);
85 PETSC_EXTERN PetscErrorCode DMForestGetComputeAdaptivitySF(DM, PetscBool *);
86 
87 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySF(DM, PetscSF *, PetscSF *);
88 
89 /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh,
90  *                                                        2 indicates typical 2:1,
91  */
92 PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt);
93 PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *);
94 
95 /* weights for repartitioning */
96 PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode);
97 PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]);
98 
99 /* weight multiplier for refinement level: useful for sub-cycling time steps */
100 PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal);
101 PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *);
102 
103 /* this process's capacity when redistributing the cells */
104 PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal);
105 PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *);
106 
107 PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM);
108 PETSC_EXTERN PetscErrorCode DMForestSetUp(DM);
109 
110 PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *);
111 PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *);
112 
113 /* miscellaneous */
114 PETSC_EXTERN PetscErrorCode DMForestTemplate(DM,MPI_Comm,DM*);
115 
116 /* type management */
117 PETSC_EXTERN PetscErrorCode DMForestRegisterType(DMType);
118 PETSC_EXTERN PetscErrorCode DMIsForest(DM,PetscBool*);
119 
120 /* p4est */
121 PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM,PetscBool *);
122 PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM,PetscBool);
123 PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM,PetscBool *);
124 PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM,PetscBool);
125 
126 #endif
127