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