xref: /petsc/include/petscdmforest.h (revision c7eeac065615148217314da71c67287c73c19b64)
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 
18c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology);
19c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *);
20c9aa14a7SToby Isaac 
21c9aa14a7SToby Isaac /* this is the coarsest possible forest: can be any DM which we can
22c9aa14a7SToby Isaac  * convert to a DMForest */
23c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM);
24c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *);
25c9aa14a7SToby Isaac 
26c9aa14a7SToby Isaac /* these are forests from which we can adapt */
27c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCoarseForest(DM, DM);
28c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCoarseForest(DM, DM *);
29c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetFineForest(DM, DM);
30c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetFineForest(DM, DM *);
31c9aa14a7SToby Isaac 
32c9aa14a7SToby Isaac /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */
33c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt);
34c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *);
35c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt);
36c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *);
37c9aa14a7SToby Isaac 
38c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt);
39c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *);
40c9aa14a7SToby Isaac 
41c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt);
42c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *);
43c9aa14a7SToby Isaac 
44c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt);
45c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *);
46c9aa14a7SToby Isaac 
473616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *);
48c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *);
49c9aa14a7SToby Isaac 
50c9aa14a7SToby Isaac /* reserve some adaptivity types */
51c9aa14a7SToby Isaac enum {DM_FOREST_KEEP = 0,
52c9aa14a7SToby Isaac       DM_FOREST_REFINE,
53c9aa14a7SToby Isaac       DM_FOREST_COARSEN,
54c9aa14a7SToby Isaac       DM_FOREST_RESERVED_ADAPTIVITY_COUNT};
55c9aa14a7SToby Isaac 
56c9aa14a7SToby Isaac /* flag each cell with an adaptivity count: should match the cell section */
57c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellAdaptivityMarkers(DM, PetscInt[], PetscCopyMode);
58*c7eeac06SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellAdaptivityMarkers(DM, PetscInt *[]);
59c9aa14a7SToby Isaac 
60c9aa14a7SToby Isaac /*J
61c9aa14a7SToby Isaac     DMForestAdaptivityStrategy - String with the name of a PETSc DMForest adaptivity strategy
62c9aa14a7SToby Isaac 
63c9aa14a7SToby Isaac    Level: intermediate
64c9aa14a7SToby Isaac 
65c9aa14a7SToby Isaac .seealso: DMForestSetType(), DMForest
66c9aa14a7SToby Isaac J*/
67c9aa14a7SToby Isaac typedef const char* DMForestAdaptivityStrategy;
68c9aa14a7SToby Isaac #define DMFORESTADAPTALL "all"
69c9aa14a7SToby Isaac #define DMFORESTADAPTANY "any"
70c9aa14a7SToby Isaac 
71c9aa14a7SToby Isaac /* how to combine: -flags         from multiple processes,
72c9aa14a7SToby Isaac  *                 -coarsen flags from multiple children
73c9aa14a7SToby Isaac  */
743616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy);
753616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *);
76c9aa14a7SToby Isaac 
77c9aa14a7SToby Isaac /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh,
78c9aa14a7SToby Isaac  *                                                        2 indicates typical 2:1,
79c9aa14a7SToby Isaac  */
80c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt);
81c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *);
82c9aa14a7SToby Isaac 
83c9aa14a7SToby Isaac /* weights for repartitioning */
84c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode);
85*c7eeac06SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]);
86c9aa14a7SToby Isaac 
87c9aa14a7SToby Isaac /* weight multiplier for refinement level: useful for sub-cycling time steps */
88c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal);
89c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *);
90c9aa14a7SToby Isaac 
91c9aa14a7SToby Isaac /* this process's capacity when redistributing the cells */
923616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal);
933616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *);
94c9aa14a7SToby Isaac 
95c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM);
96c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetUp(DM);
97c9aa14a7SToby Isaac 
98c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *);
99c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *);
100c9aa14a7SToby Isaac 
101c9aa14a7SToby Isaac #endif
102