xref: /petsc/include/petscdmforest.h (revision e2d3c13c32f46411cebb078e95a6d8080fe80715)
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 
40 /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */
41 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt);
42 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *);
43 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt);
44 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *);
45 
46 PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt);
47 PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *);
48 
49 PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt);
50 PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *);
51 
52 PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt);
53 PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *);
54 
55 PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt);
56 PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *);
57 
58 PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *);
59 PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *);
60 
61 
62 /* flag each cell with an adaptivity count: should match the cell section */
63 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityLabel(DM, const char *);
64 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityLabel(DM, const char **);
65 
66 /*J
67     DMForestAdaptivityStrategy - String with the name of a PETSc DMForest adaptivity strategy
68 
69    Level: intermediate
70 
71 .seealso: DMForestSetType(), DMForest
72 J*/
73 typedef const char* DMForestAdaptivityStrategy;
74 #define DMFORESTADAPTALL "all"
75 #define DMFORESTADAPTANY "any"
76 
77 /* how to combine: -flags         from multiple processes,
78  *                 -coarsen flags from multiple children
79  */
80 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy);
81 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *);
82 
83 PETSC_EXTERN PetscErrorCode DMForestSetComputeAdaptivitySF(DM, PetscBool);
84 PETSC_EXTERN PetscErrorCode DMForestGetComputeAdaptivitySF(DM, PetscBool *);
85 
86 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySF(DM, PetscSF *, PetscSF *);
87 
88 /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh,
89  *                                                        2 indicates typical 2:1,
90  */
91 PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt);
92 PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *);
93 
94 /* weights for repartitioning */
95 PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode);
96 PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]);
97 
98 /* weight multiplier for refinement level: useful for sub-cycling time steps */
99 PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal);
100 PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *);
101 
102 /* this process's capacity when redistributing the cells */
103 PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal);
104 PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *);
105 
106 PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM);
107 PETSC_EXTERN PetscErrorCode DMForestSetUp(DM);
108 
109 PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *);
110 PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *);
111 
112 /* miscellaneous */
113 PETSC_EXTERN PetscErrorCode DMForestTemplate(DM,MPI_Comm,DM*);
114 
115 /* type management */
116 PETSC_EXTERN PetscErrorCode DMForestRegisterType(DMType);
117 PETSC_EXTERN PetscErrorCode DMIsForest(DM,PetscBool*);
118 
119 /* p4est */
120 PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM,PetscBool *);
121 PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM,PetscBool);
122 PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM,PetscBool *);
123 PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM,PetscBool);
124 
125 #endif
126