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