1 /* 2 DMFOREST, for parallel, hierarchically refined, distributed mesh problems. 3 */ 4 #if !defined(PETSCDMFOREST_H) 5 #define PETSCDMFOREST_H 6 7 #include <petscdm.h> 8 9 /*J 10 DMForestTopology - String with the name of a PETSc DMFOREST base mesh topology. The topology is a string (e.g. 11 "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during 12 DMSetUp(). 13 14 Level: beginner 15 16 .seealso: DMForestSetTopology(), DMForestGetTopology(), DMFOREST 17 J*/ 18 typedef const char* DMForestTopology; 19 20 /* Just a name for the shape of the domain */ 21 PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology); 22 PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *); 23 24 /* this is the coarsest possible forest: can be any DM which we can 25 * convert to a DMForest (right now: plex) */ 26 PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM); 27 PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *); 28 PETSC_EXTERN PetscErrorCode DMForestSetBaseCoordinateMapping(DM, PetscErrorCode(*)(DM,PetscInt,PetscInt,const PetscReal[],PetscReal[],void*),void*); 29 PETSC_EXTERN PetscErrorCode DMForestGetBaseCoordinateMapping(DM, PetscErrorCode(**)(DM,PetscInt,PetscInt,const PetscReal[],PetscReal[],void*),void*); 30 31 /* this is the forest from which we adapt */ 32 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityForest(DM, DM); 33 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityForest(DM, DM *); 34 35 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityPurpose(DM, DMAdaptFlag); 36 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityPurpose(DM, DMAdaptFlag*); 37 38 /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */ 39 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt); 40 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *); 41 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt); 42 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *); 43 44 PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt); 45 PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *); 46 47 PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt); 48 PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *); 49 50 PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt); 51 PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *); 52 53 PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt); 54 PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *); 55 56 PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *); 57 PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *); 58 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: DMForestSetType(), DMForest 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 PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *); 110 PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *); 111 112 /* miscellaneous */ 113 PETSC_EXTERN PetscErrorCode DMForestTemplate(DM,MPI_Comm,DM*); 114 115 /* type management */ 116 PETSC_EXTERN PetscErrorCode DMForestRegisterType(DMType); 117 PETSC_EXTERN PetscErrorCode DMIsForest(DM,PetscBool*); 118 119 /* p4est */ 120 PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM,PetscBool *); 121 PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM,PetscBool); 122 PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM,PetscBool *); 123 PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM,PetscBool); 124 125 #endif 126