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 11 12 Level: beginner 13 14 .seealso: DMForestSetType(), DMForest 15 J*/ 16 typedef const char* DMForestTopology; 17 18 /* Just a name for the shape of the domain */ 19 PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology); 20 PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *); 21 22 /* this is the coarsest possible forest: can be any DM which we can 23 * convert to a DMForest (right now: plex) */ 24 PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM); 25 PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *); 26 PETSC_EXTERN PetscErrorCode DMForestSetBaseCoordinateMapping(DM, PetscErrorCode(*)(DM,PetscInt,PetscInt,const PetscReal[],PetscReal[],void*),void*); 27 PETSC_EXTERN PetscErrorCode DMForestGetBaseCoordinateMapping(DM, PetscErrorCode(**)(DM,PetscInt,PetscInt,const PetscReal[],PetscReal[],void*),void*); 28 29 /* this is the forest from which we adapt */ 30 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityForest(DM, DM); 31 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityForest(DM, DM *); 32 33 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityPurpose(DM, DMAdaptFlag); 34 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityPurpose(DM, DMAdaptFlag*); 35 36 /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */ 37 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt); 38 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *); 39 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt); 40 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *); 41 42 PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt); 43 PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *); 44 45 PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt); 46 PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *); 47 48 PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt); 49 PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *); 50 51 PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt); 52 PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *); 53 54 PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *); 55 PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *); 56 57 58 /* flag each cell with an adaptivity count: should match the cell section */ 59 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityLabel(DM, DMLabel); 60 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityLabel(DM, DMLabel*); 61 62 /*J 63 DMForestAdaptivityStrategy - String with the name of a PETSc DMForest adaptivity strategy 64 65 Level: intermediate 66 67 .seealso: DMForestSetType(), DMForest 68 J*/ 69 typedef const char* DMForestAdaptivityStrategy; 70 #define DMFORESTADAPTALL "all" 71 #define DMFORESTADAPTANY "any" 72 73 /* how to combine: -flags from multiple processes, 74 * -coarsen flags from multiple children 75 */ 76 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy); 77 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *); 78 79 PETSC_EXTERN PetscErrorCode DMForestSetComputeAdaptivitySF(DM, PetscBool); 80 PETSC_EXTERN PetscErrorCode DMForestGetComputeAdaptivitySF(DM, PetscBool *); 81 82 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySF(DM, PetscSF *, PetscSF *); 83 84 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivitySuccess(DM, PetscBool *); 85 86 PETSC_EXTERN PetscErrorCode DMForestTransferVec(DM, Vec, DM, Vec, PetscBool, PetscReal); 87 88 /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh, 89 * 2 indicates typical 2:1, 90 */ 91 PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt); 92 PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *); 93 94 /* weights for repartitioning */ 95 PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode); 96 PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]); 97 98 /* weight multiplier for refinement level: useful for sub-cycling time steps */ 99 PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal); 100 PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *); 101 102 /* this process's capacity when redistributing the cells */ 103 PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal); 104 PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *); 105 106 PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM); 107 PETSC_EXTERN PetscErrorCode DMForestSetUp(DM); 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