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