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