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