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