1c9aa14a7SToby Isaac /* 268d54884SBarry Smith DMFOREST, for parallel, hierarchically refined, distributed mesh problems. 3c9aa14a7SToby Isaac */ 4a4963045SJacob Faibussowitsch #pragma once 5c9aa14a7SToby Isaac 6c9aa14a7SToby Isaac #include <petscdm.h> 7c9aa14a7SToby Isaac 8ac09b921SBarry Smith /* SUBMANSEC = DMForest */ 9ac09b921SBarry Smith 10c9aa14a7SToby Isaac /*J 115cb80ecdSBarry Smith DMForestTopology - String with the name of a PETSc `DMFOREST` base mesh topology. The topology is a string (e.g. 125cb80ecdSBarry Smith "cube", "shell") and can be interpreted by subtypes of `DMFOREST`) to construct the base `DM` of a forest during 135cb80ecdSBarry Smith `DMSetUp()`. 14c9aa14a7SToby Isaac 15c9aa14a7SToby Isaac Level: beginner 16c9aa14a7SToby Isaac 17*af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMForestSetTopology()`, `DMForestGetTopology()`, `DMFOREST` 18c9aa14a7SToby Isaac J*/ 19c9aa14a7SToby Isaac typedef const char *DMForestTopology; 20c9aa14a7SToby Isaac 21bb2ed840SToby Isaac /* Just a name for the shape of the domain */ 22c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology); 23c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *); 24c9aa14a7SToby Isaac 25c9aa14a7SToby Isaac /* this is the coarsest possible forest: can be any DM which we can 26bb2ed840SToby Isaac * convert to a DMForest (right now: plex) */ 27c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM); 28c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *); 2999478f86SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetBaseCoordinateMapping(DM, PetscErrorCode (*)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *); 3099478f86SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetBaseCoordinateMapping(DM, PetscErrorCode (**)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *); 31c9aa14a7SToby Isaac 32ba936b91SToby Isaac /* this is the forest from which we adapt */ 33ba936b91SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityForest(DM, DM); 34ba936b91SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityForest(DM, DM *); 35c9aa14a7SToby Isaac 36a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityPurpose(DM, DMAdaptFlag); 37a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityPurpose(DM, DMAdaptFlag *); 3826d9498aSToby Isaac 39c9aa14a7SToby Isaac /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */ 40c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt); 41c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *); 42c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt); 43c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *); 44c9aa14a7SToby Isaac 45c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt); 46c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *); 47c9aa14a7SToby Isaac 48c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt); 49c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *); 50c9aa14a7SToby Isaac 51c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt); 52c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *); 53c9aa14a7SToby Isaac 5456ba9f64SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt); 5556ba9f64SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *); 5656ba9f64SToby Isaac 573616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *); 58c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *); 59c9aa14a7SToby Isaac 60c9aa14a7SToby Isaac /* flag each cell with an adaptivity count: should match the cell section */ 61a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityLabel(DM, DMLabel); 62a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityLabel(DM, DMLabel *); 63c9aa14a7SToby Isaac 64c9aa14a7SToby Isaac /*J 6587497f52SBarry Smith DMForestAdaptivityStrategy - String with the name of a PETSc `DMFOREST` adaptivity strategy 66c9aa14a7SToby Isaac 67c9aa14a7SToby Isaac Level: intermediate 68c9aa14a7SToby Isaac 69*af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMFOREST`, `DMForestSetAdaptivityStrategy()`, `DMForestGetAdaptivityStrategy()`, `DMForestSetGradeFactor()` 70c9aa14a7SToby Isaac J*/ 71c9aa14a7SToby Isaac typedef const char *DMForestAdaptivityStrategy; 72c9aa14a7SToby Isaac #define DMFORESTADAPTALL "all" 73c9aa14a7SToby Isaac #define DMFORESTADAPTANY "any" 74c9aa14a7SToby Isaac 75c9aa14a7SToby Isaac /* how to combine: -flags from multiple processes, 76c9aa14a7SToby Isaac * -coarsen flags from multiple children 77c9aa14a7SToby Isaac */ 783616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy); 793616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *); 80c9aa14a7SToby Isaac 81bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetComputeAdaptivitySF(DM, PetscBool); 82bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetComputeAdaptivitySF(DM, PetscBool *); 83bf9b5d84SToby Isaac 84bf9b5d84SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySF(DM, PetscSF *, PetscSF *); 85bf9b5d84SToby Isaac 862a133e43SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySuccess(DM, PetscBool *); 872a133e43SToby Isaac 880eb7e1eaSToby Isaac PETSC_EXTERN PetscErrorCode DMForestTransferVec(DM, Vec, DM, Vec, PetscBool, PetscReal); 89ac34a06fSStefano Zampini PETSC_EXTERN PetscErrorCode DMForestTransferVecFromBase(DM, Vec, Vec); 9080b27e07SToby Isaac 91c9aa14a7SToby Isaac /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh, 92c9aa14a7SToby Isaac * 2 indicates typical 2:1, 93c9aa14a7SToby Isaac */ 94c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt); 95c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *); 96c9aa14a7SToby Isaac 97c9aa14a7SToby Isaac /* weights for repartitioning */ 98c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode); 99c7eeac06SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]); 100c9aa14a7SToby Isaac 101c9aa14a7SToby Isaac /* weight multiplier for refinement level: useful for sub-cycling time steps */ 102c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal); 103c9aa14a7SToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *); 104c9aa14a7SToby Isaac 105c9aa14a7SToby Isaac /* this process's capacity when redistributing the cells */ 1063616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal); 1073616b93eSToby Isaac PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *); 108c9aa14a7SToby Isaac 109a0452a8eSToby Isaac /* miscellaneous */ 11020e8089bSToby Isaac PETSC_EXTERN PetscErrorCode DMForestTemplate(DM, MPI_Comm, DM *); 111a0452a8eSToby Isaac 11227d4645fSToby Isaac /* type management */ 11327d4645fSToby Isaac PETSC_EXTERN PetscErrorCode DMForestRegisterType(DMType); 11427d4645fSToby Isaac PETSC_EXTERN PetscErrorCode DMIsForest(DM, PetscBool *); 11527d4645fSToby Isaac 116a0452a8eSToby Isaac /* p4est */ 117a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM, PetscBool *); 118a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM, PetscBool); 119a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM, PetscBool *); 120a0452a8eSToby Isaac PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM, PetscBool); 121