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