1 /* 2 DMPlex, for parallel unstructured distributed mesh problems. 3 */ 4 #if !defined(PETSCDMPLEX_H) 5 #define PETSCDMPLEX_H 6 7 #include <petscsection.h> 8 #include <petscpartitioner.h> 9 #include <petscdm.h> 10 #include <petscdmplextypes.h> 11 #include <petscdt.h> 12 #include <petscfe.h> 13 #include <petscfv.h> 14 #include <petscsftypes.h> 15 #include <petscdmfield.h> 16 #include <petscviewer.h> 17 18 /* SUBMANSEC = DMPlex */ 19 20 PETSC_EXTERN PetscErrorCode PetscPartitionerDMPlexPartition(PetscPartitioner, DM, PetscSection, PetscSection, IS *); 21 22 PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellList(DM, PetscInt, PetscInt, PetscInt, const PetscInt[]); 23 PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellListParallel(DM, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscSF *, PetscInt **); 24 PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM, PetscInt, const PetscReal[]); 25 PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM, PetscInt, PetscSF, const PetscReal[]); 26 27 PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM *); 28 PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char[], PetscInt, DM *); 29 PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], DM *); 30 PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], PetscSF *, PetscInt **, DM *); 31 PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], const PetscInt[], const PetscScalar[]); 32 PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, DMPolytopeType, DM *); 33 PETSC_EXTERN PetscErrorCode DMPlexSetOptionsPrefix(DM, const char[]); 34 PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *); 35 PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt); 36 PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *); 37 PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt); 38 PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt); 39 PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]); 40 PETSC_EXTERN PetscErrorCode DMPlexGetConeTuple(DM, IS, PetscSection *, IS *); 41 PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]); 42 PETSC_EXTERN PetscErrorCode DMPlexRestoreConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]); 43 PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursiveVertices(DM, IS, IS *); 44 PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]); 45 PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt); 46 PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt); 47 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]); 48 PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]); 49 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *); 50 PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt); 51 PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]); 52 PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]); 53 PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt); 54 PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *); 55 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *); 56 PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]); 57 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 58 PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 59 PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 60 PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 61 PETSC_EXTERN PetscErrorCode DMPlexEqual(DM, DM, PetscBool *); 62 PETSC_EXTERN PetscErrorCode DMPlexOrientPoint(DM, PetscInt, PetscInt); 63 PETSC_EXTERN PetscErrorCode DMPlexOrient(DM); 64 PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 65 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM, PetscInt, PetscInt *, PetscInt *); 66 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM, PetscInt, const PetscScalar *, void *); 67 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM, PetscInt, PetscScalar *, void *); 68 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM, PetscInt, PetscInt, PetscInt *, PetscInt *); 69 PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM, PetscInt, PetscInt, PetscScalar *, void *); 70 PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM, PetscInt, PetscInt, const PetscScalar *, void *); 71 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM, PetscInt, PetscInt *, PetscInt *); 72 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM, PetscInt, const PetscScalar *, const void *); 73 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM, PetscInt, PetscScalar *, void *); 74 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM, PetscInt, PetscInt, PetscInt *, PetscInt *); 75 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM, PetscInt, PetscInt, PetscScalar *, void *); 76 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM, PetscInt, PetscInt, const PetscScalar *, void *); 77 78 /* Topological interpolation */ 79 PETSC_EXTERN const char *const DMPlexInterpolatedFlags[]; 80 81 /*E 82 DMPlexInterpolatedFlag - Describes level of topological interpolatedness. 83 84 Local or collective property depending on whether it is returned by `DMPlexIsInterpolated()` or `DMPlexIsInterpolatedCollective()`. 85 86 $ DMPLEX_INTERPOLATED_INVALID - Uninitialized value (internal use only; never returned by `DMPlexIsInterpolated()` or `DMPlexIsInterpolatedCollective()`) 87 $ DMPLEX_INTERPOLATED_NONE - Mesh is not interpolated 88 $ DMPLEX_INTERPOLATED_PARTIAL - Mesh is partially interpolated. This can e.g. mean `DMPLEX` with cells, faces and vertices but no edges represented, or a mesh with mixed cones (see `DMPlexStratify()` for an example) 89 $ DMPLEX_INTERPOLATED_MIXED - Can be returned only by `DMPlexIsInterpolatedCollective()`, meaning that `DMPlexIsInterpolated()` returns different interpolatedness on different ranks 90 $ DMPLEX_INTERPOLATED_FULL - Mesh is fully interpolated 91 92 Level: intermediate 93 94 Note: 95 An interpolated `DMPLEX` means that edges (and faces for 3d meshes) are present in the `DMPLEX` data structures. 96 97 Developer Note: 98 Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h 99 100 .seealso: `DMPLEX`, `DMPlexIsInterpolated()`, `DMPlexIsInterpolatedCollective()`, `DMPlexInterpolate()`, `DMPlexUninterpolate()` 101 E*/ 102 typedef enum { 103 DMPLEX_INTERPOLATED_INVALID = -1, 104 DMPLEX_INTERPOLATED_NONE = 0, 105 DMPLEX_INTERPOLATED_PARTIAL, 106 DMPLEX_INTERPOLATED_MIXED, 107 DMPLEX_INTERPOLATED_FULL 108 } DMPlexInterpolatedFlag; 109 110 PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 111 PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *); 112 PETSC_EXTERN PetscErrorCode DMPlexInterpolatePointSF(DM, PetscSF); 113 PETSC_EXTERN PetscErrorCode DMPlexIsInterpolated(DM, DMPlexInterpolatedFlag *); 114 PETSC_EXTERN PetscErrorCode DMPlexIsInterpolatedCollective(DM, DMPlexInterpolatedFlag *); 115 116 PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *); 117 PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 118 PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 119 PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *); 120 PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *); 121 PETSC_EXTERN PetscErrorCode DMPlexCreateLabelField(DM, DMLabel, Vec *); 122 123 PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 124 PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *); 125 PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 126 PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 127 PETSC_EXTERN PetscErrorCode DMPlexGetPointDepth(DM, PetscInt, PetscInt *); 128 PETSC_EXTERN PetscErrorCode DMPlexGetPointHeight(DM, PetscInt, PetscInt *); 129 PETSC_EXTERN PetscErrorCode DMPlexGetCellTypeLabel(DM, DMLabel *); 130 PETSC_EXTERN PetscErrorCode DMPlexGetCellType(DM, PetscInt, DMPolytopeType *); 131 PETSC_EXTERN PetscErrorCode DMPlexSetCellType(DM, PetscInt, DMPolytopeType); 132 PETSC_EXTERN PetscErrorCode DMPlexComputeCellTypes(DM); 133 PETSC_EXTERN PetscErrorCode DMPlexInvertCell(DMPolytopeType, PetscInt[]); 134 PETSC_EXTERN PetscErrorCode DMPlexReorderCell(DM, PetscInt, PetscInt[]); 135 136 /* Topological Operations */ 137 PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **); 138 PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **); 139 PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **); 140 PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **); 141 PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **); 142 PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **); 143 PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 144 PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 145 PETSC_EXTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **); 146 PETSC_EXTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **); 147 148 /*E 149 DMPlexTPSType - Type of triply-periodic surface for a `DMPLEX` 150 151 $ DMPLEX_TPS_SCHWARZ_P - Schwarz Primitive surface, defined by the equation cos(x) + cos(y) + cos(z) = 0. 152 $ DMPLEX_TPS_GYROID - Gyroid surface, defined by the equation sin(x)cos(y) + sin(y)cos(z) + sin(z)cos(x) = 0 153 154 Level: intermediate 155 156 Developer Note: 157 Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h 158 159 .seealso: `DMPlexCreateTPSMesh()` 160 E*/ 161 typedef enum { 162 DMPLEX_TPS_SCHWARZ_P, 163 DMPLEX_TPS_GYROID 164 } DMPlexTPSType; 165 PETSC_EXTERN const char *const DMPlexTPSTypes[]; 166 167 /* Mesh Generation */ 168 PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char[], PetscBool, DM *); 169 PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 170 PETSC_EXTERN PetscErrorCode DMPlexCreateCoordinateSpace(DM, PetscInt, void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 171 PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscReal, DM *); 172 PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *); 173 PETSC_EXTERN PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscBool, DM *); 174 PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, PetscReal, DM *); 175 PETSC_EXTERN PetscErrorCode DMPlexCreateBallMesh(MPI_Comm, PetscInt, PetscReal, DM *); 176 PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, DMBoundaryType, DM *); 177 PETSC_EXTERN PetscErrorCode DMPlexCreateTPSMesh(MPI_Comm, DMPlexTPSType, const PetscInt[], const DMBoundaryType[], PetscBool, PetscInt, PetscInt, PetscReal, DM *); 178 PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *); 179 PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *); 180 PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, const PetscReal[], const PetscReal[], DM *); 181 PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 182 PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModel(DM); 183 184 PETSC_EXTERN PetscErrorCode DMPlexCheck(DM); 185 PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM); 186 PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscInt); 187 PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscInt); 188 PETSC_EXTERN PetscErrorCode DMPlexCheckGeometry(DM); 189 PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM, PetscSF); 190 PETSC_EXTERN PetscErrorCode DMPlexCheckInterfaceCones(DM); 191 PETSC_EXTERN PetscErrorCode DMPlexCheckCellShape(DM, PetscBool, PetscReal); 192 PETSC_EXTERN PetscErrorCode DMPlexComputeOrthogonalQuality(DM, PetscFV, PetscReal, Vec *, DMLabel *); 193 194 PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *); 195 PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *); 196 197 PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], const char[], PetscBool, DM *); 198 PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 199 PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char[], PetscBool, DM *); 200 PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 201 PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *); 202 PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *); 203 PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char[], PetscBool, DM *); 204 PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *); 205 PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char[], PetscBool, DM *); 206 PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char[], PetscBool, DM *); 207 PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char[], PetscBool, DM *); 208 PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm, const char[], DM *); 209 PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSLiteFromFile(MPI_Comm, const char[], DM *); 210 211 PETSC_EXTERN PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo); 212 PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetId(PetscViewer, int *); 213 PETSC_EXTERN PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer, PetscInt); 214 PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer, PetscInt *); 215 216 /* Mesh Partitioning and Distribution */ 217 #define DMPLEX_OVERLAP_MANUAL -1 218 219 PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 220 PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *); 221 PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner); 222 PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *); 223 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel); 224 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel); 225 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel); 226 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel); 227 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *); 228 PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool); 229 PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *); 230 PETSC_EXTERN PetscErrorCode DMPlexIsDistributed(DM, PetscBool *); 231 PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF *, DM *); 232 PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *); 233 PETSC_EXTERN PetscErrorCode DMPlexGetOverlap(DM, PetscInt *); 234 PETSC_EXTERN PetscErrorCode DMPlexSetOverlap(DM, DM, PetscInt); 235 PETSC_EXTERN PetscErrorCode DMPlexDistributeGetDefault(DM, PetscBool *); 236 PETSC_EXTERN PetscErrorCode DMPlexDistributeSetDefault(DM, PetscBool); 237 PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM, PetscSF, PetscSection, Vec, PetscSection, Vec); 238 PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *); 239 PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM, PetscSF, PetscSection, MPI_Datatype, void *, PetscSection, void **); 240 PETSC_EXTERN PetscErrorCode DMPlexRebalanceSharedPoints(DM, PetscInt, PetscBool, PetscBool, PetscBool *); 241 PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM); 242 PETSC_EXTERN PetscErrorCode DMPlexGetGatherDM(DM, PetscSF *, DM *); 243 PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, PetscSF *, DM *); 244 PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM, PetscInt, PetscInt *, PetscInt[], void *), void *); 245 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM, PetscInt, PetscInt *, PetscInt[], void *), void **); 246 PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM, PetscBool); 247 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM, PetscBool *); 248 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]); 249 PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF); 250 PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *); 251 PETSC_EXTERN PetscErrorCode DMPlexDistributionSetName(DM, const char[]); 252 PETSC_EXTERN PetscErrorCode DMPlexDistributionGetName(DM, const char *[]); 253 254 /*E 255 DMPlexReorderDefaultFlag - Flag indicating whether the `DMPLEX` should be reordered by default 256 257 $ DMPLEX_REORDER_DEFAULT_NOTSET - Flag not set. 258 $ DMPLEX_REORDER_DEFAULT_FALSE - Do not reorder by default. 259 $ DMPLEX_REORDER_DEFAULT_TRUE - Reorder by default. 260 261 Level: intermediate 262 263 Developer Note: 264 Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h 265 266 Could be replaced with `PETSC_BOOL3` 267 268 .seealso: `DMPlexReorderSetDefault()`, `DMPlexReorderGetDefault()`, `DMPlexGetOrdering()`, `DMPlexPermute()` 269 E*/ 270 typedef enum { 271 DMPLEX_REORDER_DEFAULT_NOTSET = -1, 272 DMPLEX_REORDER_DEFAULT_FALSE = 0, 273 DMPLEX_REORDER_DEFAULT_TRUE 274 } DMPlexReorderDefaultFlag; 275 PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *); 276 PETSC_EXTERN PetscErrorCode DMPlexGetOrdering1D(DM, IS *); 277 PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *); 278 PETSC_EXTERN PetscErrorCode DMPlexReorderGetDefault(DM, DMPlexReorderDefaultFlag *); 279 PETSC_EXTERN PetscErrorCode DMPlexReorderSetDefault(DM, DMPlexReorderDefaultFlag); 280 281 PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *); 282 PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *); 283 PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *); 284 PETSC_EXTERN PetscErrorCode DMPlexCreatePointSF(DM, PetscSF, PetscBool, PetscSF *); 285 PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabel(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *); 286 PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabelFromLabels(DM, PetscInt, const DMLabel[], const PetscInt[], PetscInt, const DMLabel[], const PetscInt[], PetscSection, IS, PetscSection, IS, DMLabel *); 287 PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *); 288 PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *); 289 290 /* Submesh Support */ 291 PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM *); 292 PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, PetscInt, DMLabel *, DMLabel *, DM *, DM *); 293 PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel *); 294 PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 295 PETSC_EXTERN PetscErrorCode DMPlexGetSubpointIS(DM, IS *); 296 297 PETSC_EXTERN PetscErrorCode DMGetEnclosureRelation(DM, DM, DMEnclosureType *); 298 PETSC_EXTERN PetscErrorCode DMGetEnclosurePoint(DM, DM, DMEnclosureType, PetscInt, PetscInt *); 299 300 PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 301 PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscInt, PetscBool, DM); 302 PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel); 303 PETSC_EXTERN PetscErrorCode DMPlexLabelAddFaceCells(DM, DMLabel); 304 PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel); 305 306 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 307 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 308 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 309 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 310 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal[], PetscReal *)); 311 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal[], PetscReal *)); 312 PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *); 313 PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *); 314 PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool); 315 316 /* Support for cell-vertex meshes */ 317 PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *); 318 PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt[], PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscBool *); 319 320 /* Geometry Support */ 321 PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *); 322 PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal); 323 PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]); 324 PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]); 325 PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]); 326 327 /* Point Location */ 328 typedef struct _PetscGridHash *PetscGridHash; 329 PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *); 330 PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar[]); 331 PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt[], const PetscReal[]); 332 PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar[], PetscInt[], PetscInt[]); 333 PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *); 334 PETSC_EXTERN PetscErrorCode DMPlexFindVertices(DM, Vec, PetscReal, IS *); 335 336 /* FVM Support */ 337 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal[], PetscReal[]); 338 PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *); 339 PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *); 340 PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *); 341 342 /* FEM Support */ 343 PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *); 344 PETSC_EXTERN PetscErrorCode DMPlexGetGeometryFVM(DM, Vec *, Vec *, PetscReal *); 345 PETSC_EXTERN PetscErrorCode DMPlexGetGradientDM(DM, PetscFV, DM *); 346 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec); 347 PETSC_EXTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec); 348 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec); 349 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void *, Vec); 350 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void *, Vec); 351 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM, PetscReal, Vec, Vec, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], PetscErrorCode (*)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *), void *, Vec); 352 PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel); 353 354 PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, DMLabel[], const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const IS[], const IS[], IS, PetscSection *); 355 PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection *); 356 357 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 358 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 359 PETSC_EXTERN PetscErrorCode DMPlexGetCellCoordinates(DM, PetscInt, PetscBool *, PetscInt *, const PetscScalar *[], PetscScalar *[]); 360 PETSC_EXTERN PetscErrorCode DMPlexRestoreCellCoordinates(DM, PetscInt, PetscBool *, PetscInt *, const PetscScalar *[], PetscScalar *[]); 361 PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[]); 362 PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[]); 363 PETSC_EXTERN PetscErrorCode DMPlexShearGeometry(DM, DMDirection, PetscReal[]); 364 PETSC_EXTERN PetscErrorCode DMPlexRemapGeometry(DM, PetscReal, void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 365 366 PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 367 PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 368 PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 369 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 370 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureGeneral(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 371 PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]); 372 PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]); 373 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 374 PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]); 375 PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection); 376 PETSC_EXTERN PetscErrorCode DMPlexSetClosurePermutationTensor(DM, PetscInt, PetscSection); 377 378 PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char[], PetscInt *, DM *); 379 PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *); 380 381 PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 382 PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 383 PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 384 PETSC_EXTERN PetscErrorCode DMPlexGetGhostCellStratum(DM, PetscInt *, PetscInt *); 385 PETSC_EXTERN PetscErrorCode DMPlexGetSimplexOrBoxCells(DM, PetscInt, PetscInt *, PetscInt *); 386 PETSC_EXTERN PetscErrorCode DMPlexIsSimplex(DM, PetscBool *); 387 388 PETSC_EXTERN PetscErrorCode DMPlexGetCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **); 389 PETSC_EXTERN PetscErrorCode DMPlexRestoreCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **); 390 PETSC_EXTERN PetscErrorCode DMPlexGetFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **); 391 PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **); 392 PETSC_EXTERN PetscErrorCode DMPlexGetFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **); 393 PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **); 394 395 PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 396 PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 397 398 typedef struct { 399 DM dm; 400 Vec u; /* The base vector for the Jacobian action J(u) x */ 401 Mat J; /* Preconditioner for testing */ 402 void *user; 403 } JacActionCtx; 404 405 PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt); 406 PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt *); 407 PETSC_EXTERN PetscErrorCode DMPlexGetActivePoint(DM, PetscInt *); 408 PETSC_EXTERN PetscErrorCode DMPlexSetActivePoint(DM, PetscInt); 409 PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscScalar[]), InsertMode, Vec); 410 PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 411 PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]); 412 PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, Vec); 413 PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *); 414 PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *); 415 PETSC_EXTERN PetscErrorCode DMPlexComputeBdIntegral(DM, Vec, DMLabel, PetscInt, const PetscInt[], void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscScalar *, void *); 416 PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, PetscBool, Mat, void *); 417 PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *); 418 PETSC_EXTERN PetscErrorCode DMPlexComputeClementInterpolant(DM, Vec, Vec); 419 PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec); 420 PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *); 421 PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *); 422 PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *); 423 424 PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, PetscInt, MatNullSpace *); 425 PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *); 426 427 PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM, void *, void *, void *); 428 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *); 429 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *); 430 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *); 431 PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, PetscWeakForm, PetscFormKey, Vec, Vec, Vec); 432 PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, PetscWeakForm, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat); 433 434 PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *); 435 PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 436 PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *); 437 PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 438 PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM, PetscReal, Vec, Vec, void *); 439 440 PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 441 PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM, Vec, Vec); 442 443 /* anchors */ 444 PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection *, IS *); 445 PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS); 446 /* tree */ 447 PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM); 448 PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM *); 449 PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *); 450 PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM *); 451 PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt[], PetscInt[]); 452 PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]); 453 PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *); 454 PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]); 455 PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *); 456 PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *); 457 PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM, Vec, DM, Vec, PetscSF, PetscSF, PetscInt *, PetscInt *, PetscBool, PetscReal); 458 459 PETSC_EXTERN PetscErrorCode DMPlexMonitorThroughput(DM, void *); 460 461 /* natural order */ 462 PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *); 463 PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF); 464 PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *); 465 PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec); 466 PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec); 467 PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec); 468 PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec); 469 PETSC_EXTERN PetscErrorCode DMPlexCreateNaturalVector(DM, Vec *); 470 471 /* mesh adaptation */ 472 PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char[], DM *); 473 PETSC_EXTERN PetscErrorCode DMPlexSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]); 474 PETSC_EXTERN PetscErrorCode DMPlexMetricSetFromOptions(DM); 475 PETSC_EXTERN PetscErrorCode DMPlexMetricSetIsotropic(DM, PetscBool); 476 PETSC_EXTERN PetscErrorCode DMPlexMetricIsIsotropic(DM, PetscBool *); 477 PETSC_EXTERN PetscErrorCode DMPlexMetricSetUniform(DM, PetscBool); 478 PETSC_EXTERN PetscErrorCode DMPlexMetricIsUniform(DM, PetscBool *); 479 PETSC_EXTERN PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM, PetscBool); 480 PETSC_EXTERN PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM, PetscBool *); 481 PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoInsertion(DM, PetscBool); 482 PETSC_EXTERN PetscErrorCode DMPlexMetricNoInsertion(DM, PetscBool *); 483 PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSwapping(DM, PetscBool); 484 PETSC_EXTERN PetscErrorCode DMPlexMetricNoSwapping(DM, PetscBool *); 485 PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoMovement(DM, PetscBool); 486 PETSC_EXTERN PetscErrorCode DMPlexMetricNoMovement(DM, PetscBool *); 487 PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSurf(DM, PetscBool); 488 PETSC_EXTERN PetscErrorCode DMPlexMetricNoSurf(DM, PetscBool *); 489 PETSC_EXTERN PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM, PetscReal); 490 PETSC_EXTERN PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM, PetscReal *); 491 PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM, PetscReal); 492 PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM, PetscReal *); 493 PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM, PetscReal); 494 PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM, PetscReal *); 495 PETSC_EXTERN PetscErrorCode DMPlexMetricSetTargetComplexity(DM, PetscReal); 496 PETSC_EXTERN PetscErrorCode DMPlexMetricGetTargetComplexity(DM, PetscReal *); 497 PETSC_EXTERN PetscErrorCode DMPlexMetricSetNormalizationOrder(DM, PetscReal); 498 PETSC_EXTERN PetscErrorCode DMPlexMetricGetNormalizationOrder(DM, PetscReal *); 499 PETSC_EXTERN PetscErrorCode DMPlexMetricSetGradationFactor(DM, PetscReal); 500 PETSC_EXTERN PetscErrorCode DMPlexMetricGetGradationFactor(DM, PetscReal *); 501 PETSC_EXTERN PetscErrorCode DMPlexMetricSetHausdorffNumber(DM, PetscReal); 502 PETSC_EXTERN PetscErrorCode DMPlexMetricGetHausdorffNumber(DM, PetscReal *); 503 PETSC_EXTERN PetscErrorCode DMPlexMetricSetVerbosity(DM, PetscInt); 504 PETSC_EXTERN PetscErrorCode DMPlexMetricGetVerbosity(DM, PetscInt *); 505 PETSC_EXTERN PetscErrorCode DMPlexMetricSetNumIterations(DM, PetscInt); 506 PETSC_EXTERN PetscErrorCode DMPlexMetricGetNumIterations(DM, PetscInt *); 507 PETSC_EXTERN PetscErrorCode DMPlexMetricCreate(DM, PetscInt, Vec *); 508 PETSC_EXTERN PetscErrorCode DMPlexMetricCreateUniform(DM, PetscInt, PetscReal, Vec *); 509 PETSC_EXTERN PetscErrorCode DMPlexMetricCreateIsotropic(DM, PetscInt, Vec, Vec *); 510 PETSC_EXTERN PetscErrorCode DMPlexMetricDeterminantCreate(DM, PetscInt, Vec *, DM *); 511 PETSC_EXTERN PetscErrorCode DMPlexMetricEnforceSPD(DM, Vec, PetscBool, PetscBool, Vec, Vec); 512 PETSC_EXTERN PetscErrorCode DMPlexMetricNormalize(DM, Vec, PetscBool, PetscBool, Vec, Vec); 513 PETSC_EXTERN PetscErrorCode DMPlexMetricAverage(DM, PetscInt, PetscReal[], Vec[], Vec); 514 PETSC_EXTERN PetscErrorCode DMPlexMetricAverage2(DM, Vec, Vec, Vec); 515 PETSC_EXTERN PetscErrorCode DMPlexMetricAverage3(DM, Vec, Vec, Vec, Vec); 516 PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection(DM, PetscInt, Vec[], Vec); 517 PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection2(DM, Vec, Vec, Vec); 518 PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection3(DM, Vec, Vec, Vec, Vec); 519 520 PETSC_EXTERN PetscErrorCode DMPlexGlobalToLocalBasis(DM, Vec); 521 PETSC_EXTERN PetscErrorCode DMPlexLocalToGlobalBasis(DM, Vec); 522 PETSC_EXTERN PetscErrorCode DMPlexCreateBasisRotation(DM, PetscReal, PetscReal, PetscReal); 523 524 /* storage version */ 525 #define DMPLEX_STORAGE_VERSION_STABLE "1.0.0" 526 #define DMPLEX_STORAGE_VERSION_LATEST "2.0.0" 527 528 PETSC_EXTERN PetscErrorCode DMPlexTopologyView(DM, PetscViewer); 529 PETSC_EXTERN PetscErrorCode DMPlexCoordinatesView(DM, PetscViewer); 530 PETSC_EXTERN PetscErrorCode DMPlexLabelsView(DM, PetscViewer); 531 PETSC_EXTERN PetscErrorCode DMPlexSectionView(DM, PetscViewer, DM); 532 PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorView(DM, PetscViewer, DM, Vec); 533 PETSC_EXTERN PetscErrorCode DMPlexLocalVectorView(DM, PetscViewer, DM, Vec); 534 PETSC_EXTERN PetscErrorCode DMPlexVecView1D(DM, PetscInt, Vec[], PetscViewer); 535 536 PETSC_EXTERN PetscErrorCode DMPlexTopologyLoad(DM, PetscViewer, PetscSF *); 537 PETSC_EXTERN PetscErrorCode DMPlexCoordinatesLoad(DM, PetscViewer, PetscSF); 538 PETSC_EXTERN PetscErrorCode DMPlexLabelsLoad(DM, PetscViewer, PetscSF); 539 PETSC_EXTERN PetscErrorCode DMPlexSectionLoad(DM, PetscViewer, DM, PetscSF, PetscSF *, PetscSF *); 540 PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorLoad(DM, PetscViewer, DM, PetscSF, Vec); 541 PETSC_EXTERN PetscErrorCode DMPlexLocalVectorLoad(DM, PetscViewer, DM, PetscSF, Vec); 542 543 PETSC_EXTERN PetscErrorCode DMPlexGetLocalOffsets(DM, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt **); 544 545 /* point queue */ 546 PETSC_EXTERN PetscErrorCode DMPlexPointQueueCreate(PetscInt, DMPlexPointQueue *); 547 PETSC_EXTERN PetscErrorCode DMPlexPointQueueDestroy(DMPlexPointQueue *); 548 PETSC_EXTERN PetscErrorCode DMPlexPointQueueEnsureSize(DMPlexPointQueue); 549 PETSC_EXTERN PetscErrorCode DMPlexPointQueueEnqueue(DMPlexPointQueue, PetscInt); 550 PETSC_EXTERN PetscErrorCode DMPlexPointQueueDequeue(DMPlexPointQueue, PetscInt *); 551 PETSC_EXTERN PetscErrorCode DMPlexPointQueueFront(DMPlexPointQueue, PetscInt *); 552 PETSC_EXTERN PetscErrorCode DMPlexPointQueueBack(DMPlexPointQueue, PetscInt *); 553 PETSC_EXTERN PetscBool DMPlexPointQueueEmpty(DMPlexPointQueue); 554 PETSC_EXTERN PetscErrorCode DMPlexPointQueueEmptyCollective(PetscObject, DMPlexPointQueue, PetscBool *); 555 556 #endif 557