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