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 3126d9498aSToby Isaac /* reserve some adaptivity types */ 3226d9498aSToby Isaac enum {DM_FOREST_KEEP = 0, 3326d9498aSToby Isaac DM_FOREST_REFINE, 3426d9498aSToby Isaac DM_FOREST_COARSEN, 3526d9498aSToby Isaac DM_FOREST_RESERVED_ADAPTIVITY_COUNT}; 3626d9498aSToby Isaac 3726d9498aSToby Isaac typedef PetscInt DMForestAdaptivityPurpose; 3826d9498aSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityPurpose(DM, DMForestAdaptivityPurpose); 3926d9498aSToby 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 83*bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetComputeAdaptivitySF(DM, PetscBool); 84*bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetComputeAdaptivitySF(DM, PetscBool *); 85*bf9b5d84SToby Isaac 86*bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySF(DM, PetscSF *, PetscSF *); 87*bf9b5d84SToby Isaac 88c9aa14a7SToby Isaac /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh, 89c9aa14a7SToby Isaac * 2 indicates typical 2:1, 90c9aa14a7SToby Isaac */ 91c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt); 92c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *); 93c9aa14a7SToby Isaac 94c9aa14a7SToby Isaac /* weights for repartitioning */ 95c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode); 96c7eeac06SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]); 97c9aa14a7SToby Isaac 98c9aa14a7SToby Isaac /* weight multiplier for refinement level: useful for sub-cycling time steps */ 99c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal); 100c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *); 101c9aa14a7SToby Isaac 102c9aa14a7SToby Isaac /* this process's capacity when redistributing the cells */ 1033616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal); 1043616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *); 105c9aa14a7SToby Isaac 106c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM); 107c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetUp(DM); 108c9aa14a7SToby Isaac 109c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *); 110c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *); 111c9aa14a7SToby Isaac 112a0452a8eSToby Isaac /* miscellaneous */ 11320e8089bSToby Isaac PETSC_EXTERN PetscErrorCode DMForestTemplate(DM,MPI_Comm,DM*); 114a0452a8eSToby Isaac 11527d4645fSToby Isaac /* type management */ 11627d4645fSToby Isaac PETSC_EXTERN PetscErrorCode DMForestRegisterType(DMType); 11727d4645fSToby Isaac PETSC_EXTERN PetscErrorCode DMIsForest(DM,PetscBool*); 11827d4645fSToby Isaac 119a0452a8eSToby Isaac /* p4est */ 120a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM,PetscBool *); 121a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM,PetscBool); 122a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM,PetscBool *); 123a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM,PetscBool); 124a0452a8eSToby Isaac 125c9aa14a7SToby Isaac #endif 126