xref: /petsc/include/petscdmforest.h (revision 0eb7e1eae252ea4368f51905e6ed2ea3dd84b0fd)
1c9aa14a7SToby Isaac /*
2c9aa14a7SToby Isaac   DMForest, for parallel, hierarchically refined, distributed mesh problems.
3c9aa14a7SToby Isaac */
4c9aa14a7SToby Isaac #if !defined(__PETSCDMFOREST_H)
5c9aa14a7SToby Isaac #define __PETSCDMFOREST_H
6c9aa14a7SToby Isaac 
7c9aa14a7SToby Isaac #include <petscdm.h>
8c9aa14a7SToby Isaac 
9c9aa14a7SToby Isaac /*J
10c9aa14a7SToby Isaac     DMForestTopology - String with the name of a PETSc DMForest base mesh topology
11c9aa14a7SToby Isaac 
12c9aa14a7SToby Isaac    Level: beginner
13c9aa14a7SToby Isaac 
14c9aa14a7SToby Isaac .seealso: DMForestSetType(), DMForest
15c9aa14a7SToby Isaac J*/
16c9aa14a7SToby Isaac typedef const char* DMForestTopology;
17c9aa14a7SToby Isaac 
18bb2ed840SToby Isaac /* Just a name for the shape of the domain */
19c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology);
20c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *);
21c9aa14a7SToby Isaac 
22c9aa14a7SToby Isaac /* this is the coarsest possible forest: can be any DM which we can
23bb2ed840SToby Isaac  * convert to a DMForest (right now: plex) */
24c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM);
25c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *);
2699478f86SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetBaseCoordinateMapping(DM, PetscErrorCode(*)(DM,PetscInt,PetscInt,const PetscReal[],PetscReal[],void*),void*);
2799478f86SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetBaseCoordinateMapping(DM, PetscErrorCode(**)(DM,PetscInt,PetscInt,const PetscReal[],PetscReal[],void*),void*);
28c9aa14a7SToby Isaac 
29ba936b91SToby Isaac /* this is the forest from which we adapt */
30ba936b91SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityForest(DM, DM);
31ba936b91SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityForest(DM, DM *);
32c9aa14a7SToby Isaac 
3326d9498aSToby Isaac /* reserve some adaptivity types */
3426d9498aSToby Isaac enum {DM_FOREST_KEEP = 0,
3526d9498aSToby Isaac       DM_FOREST_REFINE,
3626d9498aSToby Isaac       DM_FOREST_COARSEN,
3726d9498aSToby Isaac       DM_FOREST_RESERVED_ADAPTIVITY_COUNT};
3826d9498aSToby Isaac 
3926d9498aSToby Isaac typedef PetscInt DMForestAdaptivityPurpose;
4026d9498aSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityPurpose(DM, DMForestAdaptivityPurpose);
4156c0450aSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityPurpose(DM, DMForestAdaptivityPurpose*);
4226d9498aSToby Isaac 
43c9aa14a7SToby Isaac /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */
44c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt);
45c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *);
46c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt);
47c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *);
48c9aa14a7SToby Isaac 
49c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt);
50c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *);
51c9aa14a7SToby Isaac 
52c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt);
53c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *);
54c9aa14a7SToby Isaac 
55c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt);
56c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *);
57c9aa14a7SToby Isaac 
5856ba9f64SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt);
5956ba9f64SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *);
6056ba9f64SToby Isaac 
613616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *);
62c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *);
63c9aa14a7SToby Isaac 
64c9aa14a7SToby Isaac 
65c9aa14a7SToby Isaac /* flag each cell with an adaptivity count: should match the cell section */
66ebdf65a2SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityLabel(DM, const char *);
67ba936b91SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityLabel(DM, const char **);
68c9aa14a7SToby Isaac 
69c9aa14a7SToby Isaac /*J
70c9aa14a7SToby Isaac     DMForestAdaptivityStrategy - String with the name of a PETSc DMForest adaptivity strategy
71c9aa14a7SToby Isaac 
72c9aa14a7SToby Isaac    Level: intermediate
73c9aa14a7SToby Isaac 
74c9aa14a7SToby Isaac .seealso: DMForestSetType(), DMForest
75c9aa14a7SToby Isaac J*/
76c9aa14a7SToby Isaac typedef const char* DMForestAdaptivityStrategy;
77c9aa14a7SToby Isaac #define DMFORESTADAPTALL "all"
78c9aa14a7SToby Isaac #define DMFORESTADAPTANY "any"
79c9aa14a7SToby Isaac 
80c9aa14a7SToby Isaac /* how to combine: -flags         from multiple processes,
81c9aa14a7SToby Isaac  *                 -coarsen flags from multiple children
82c9aa14a7SToby Isaac  */
833616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy);
843616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *);
85c9aa14a7SToby Isaac 
86bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetComputeAdaptivitySF(DM, PetscBool);
87bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetComputeAdaptivitySF(DM, PetscBool *);
88bf9b5d84SToby Isaac 
89bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySF(DM, PetscSF *, PetscSF *);
90bf9b5d84SToby Isaac 
91*0eb7e1eaSToby Isaac PETSC_EXTERN PetscErrorCode DMForestTransferVec(DM, Vec, DM, Vec, PetscBool, PetscReal);
9280b27e07SToby Isaac 
93c9aa14a7SToby Isaac /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh,
94c9aa14a7SToby Isaac  *                                                        2 indicates typical 2:1,
95c9aa14a7SToby Isaac  */
96c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt);
97c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *);
98c9aa14a7SToby Isaac 
99c9aa14a7SToby Isaac /* weights for repartitioning */
100c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode);
101c7eeac06SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]);
102c9aa14a7SToby Isaac 
103c9aa14a7SToby Isaac /* weight multiplier for refinement level: useful for sub-cycling time steps */
104c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal);
105c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *);
106c9aa14a7SToby Isaac 
107c9aa14a7SToby Isaac /* this process's capacity when redistributing the cells */
1083616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal);
1093616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *);
110c9aa14a7SToby Isaac 
111c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM);
112c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetUp(DM);
113c9aa14a7SToby Isaac 
114c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *);
115c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *);
116c9aa14a7SToby Isaac 
117a0452a8eSToby Isaac /* miscellaneous */
11820e8089bSToby Isaac PETSC_EXTERN PetscErrorCode DMForestTemplate(DM,MPI_Comm,DM*);
119a0452a8eSToby Isaac 
12027d4645fSToby Isaac /* type management */
12127d4645fSToby Isaac PETSC_EXTERN PetscErrorCode DMForestRegisterType(DMType);
12227d4645fSToby Isaac PETSC_EXTERN PetscErrorCode DMIsForest(DM,PetscBool*);
12327d4645fSToby Isaac 
124a0452a8eSToby Isaac /* p4est */
125a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM,PetscBool *);
126a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM,PetscBool);
127a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM,PetscBool *);
128a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM,PetscBool);
129a0452a8eSToby Isaac 
130c9aa14a7SToby Isaac #endif
131