xref: /petsc/include/petscdmplex.h (revision 1fcf445a439848b2e538efa39b3f63deb630163e)
1552f7358SJed Brown /*
2552f7358SJed Brown   DMPlex, for parallel unstructured distributed mesh problems.
3552f7358SJed Brown */
46524c165SJacob Faibussowitsch #ifndef PETSCDMPLEX_H
526bd1501SBarry Smith #define PETSCDMPLEX_H
6552f7358SJed Brown 
7224ef0b1SMatthew Knepley #include <petscsection.h>
8abe9303eSLisandro Dalcin #include <petscpartitioner.h>
9552f7358SJed Brown #include <petscdm.h>
1086fe8405SMatthew G. Knepley #include <petscdmplextypes.h>
11a0845e3aSMatthew G. Knepley #include <petscdt.h>
12a0845e3aSMatthew G. Knepley #include <petscfe.h>
130bacfa0aSMatthew G. Knepley #include <petscfv.h>
14799f573fSMatthew G. Knepley #include <petscsftypes.h>
15c330f8ffSToby Isaac #include <petscdmfield.h>
16ea8e1828Sksagiyam #include <petscviewer.h>
17552f7358SJed Brown 
18ac09b921SBarry Smith /* SUBMANSEC = DMPlex */
19ac09b921SBarry Smith 
203c41b853SStefano Zampini PETSC_EXTERN PetscErrorCode PetscPartitionerDMPlexPartition(PetscPartitioner, DM, PetscSection, PetscSection, IS *);
2177623264SMatthew G. Knepley 
225e488331SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellList(DM, PetscInt, PetscInt, PetscInt, const PetscInt[]);
23be8c289dSNicolas Barral PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellListParallel(DM, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscSF *, PetscInt **);
241edcf0b2SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM, PetscInt, const PetscReal[]);
251edcf0b2SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM, PetscInt, PetscSF, const PetscReal[]);
26b09969d6SVaclav Hapla 
27552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM *);
2827c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char[], PetscInt, DM *);
29a4a685f2SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], DM *);
30be8c289dSNicolas Barral PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], PetscSF *, PetscInt **, DM *);
3106bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], const PetscInt[], const PetscScalar[]);
329318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, DMPolytopeType, DM *);
33a9074c1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetOptionsPrefix(DM, const char[]);
34552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *);
35552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt);
36552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *);
37552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt);
38552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
390ce7577fSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeTuple(DM, IS, PetscSection *, IS *);
40af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
41af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexRestoreConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
42af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursiveVertices(DM, IS, IS *);
439f4ada15SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedCone(DM, PetscInt, const PetscInt *[], const PetscInt *[]);
449f4ada15SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreOrientedCone(DM, PetscInt, const PetscInt *[], const PetscInt *[]);
45552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
469f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
4777c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
48552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
49552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
50552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
51552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
52552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
53552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
54552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
55552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
568cb4d582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *);
57552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
58552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
59552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
61552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
624e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexEqual(DM, DM, PetscBool *);
63b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrientPoint(DM, PetscInt, PetscInt);
64d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM);
65e101f074SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
66552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM, PetscInt, PetscInt *, PetscInt *);
67081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM, PetscInt, const PetscScalar *, void *);
68a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM, PetscInt, PetscScalar *, void *);
69a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM, PetscInt, PetscInt, PetscInt *, PetscInt *);
70a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM, PetscInt, PetscInt, PetscScalar *, void *);
71081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM, PetscInt, PetscInt, const PetscScalar *, void *);
72552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM, PetscInt, PetscInt *, PetscInt *);
73552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM, PetscInt, const PetscScalar *, const void *);
74a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM, PetscInt, PetscScalar *, void *);
75a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM, PetscInt, PetscInt, PetscInt *, PetscInt *);
76a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM, PetscInt, PetscInt, PetscScalar *, void *);
77081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM, PetscInt, PetscInt, const PetscScalar *, void *);
78552f7358SJed Brown 
79a7748839SVaclav Hapla /* Topological interpolation */
80a7748839SVaclav Hapla PETSC_EXTERN const char *const DMPlexInterpolatedFlags[];
81a1cc8097SBarry Smith 
829ade3219SVaclav Hapla /*E
839ade3219SVaclav Hapla    DMPlexInterpolatedFlag - Describes level of topological interpolatedness.
849ade3219SVaclav Hapla 
8587497f52SBarry Smith    Local or collective property depending on whether it is returned by `DMPlexIsInterpolated()` or `DMPlexIsInterpolatedCollective()`.
8687497f52SBarry Smith 
8716a05f60SBarry Smith    Values:
8816a05f60SBarry Smith +  `DMPLEX_INTERPOLATED_INVALID` - Uninitialized value (internal use only; never returned by `DMPlexIsInterpolated()` or `DMPlexIsInterpolatedCollective()`)
8916a05f60SBarry Smith .  `DMPLEX_INTERPOLATED_NONE`    - Mesh is not interpolated
9016a05f60SBarry Smith .  `DMPLEX_INTERPOLATED_PARTIAL` - Mesh is partially interpolated. This can e.g. mean `DMPLEX` with cells, faces and vertices but no edges represented,
9116a05f60SBarry Smith                                    or a mesh with mixed cones (see `DMPlexStratify()` for an example)
9216a05f60SBarry Smith .  `DMPLEX_INTERPOLATED_MIXED`   - Can be returned only by `DMPlexIsInterpolatedCollective()`, meaning that `DMPlexIsInterpolated()` returns different interpolatedness on different ranks
9316a05f60SBarry Smith -  `DMPLEX_INTERPOLATED_FULL`    - Mesh is fully interpolated
949ade3219SVaclav Hapla 
959ade3219SVaclav Hapla    Level: intermediate
969ade3219SVaclav Hapla 
9787497f52SBarry Smith    Note:
9887497f52SBarry Smith    An interpolated `DMPLEX` means that edges (and faces for 3d meshes) are present in the `DMPLEX` data structures.
9987497f52SBarry Smith 
1009ade3219SVaclav Hapla    Developer Note:
1019ade3219SVaclav Hapla    Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
1029ade3219SVaclav Hapla 
10387497f52SBarry Smith .seealso: `DMPLEX`, `DMPlexIsInterpolated()`, `DMPlexIsInterpolatedCollective()`, `DMPlexInterpolate()`, `DMPlexUninterpolate()`
1049ade3219SVaclav Hapla E*/
105a7748839SVaclav Hapla typedef enum {
1067d0f5628SVaclav Hapla   DMPLEX_INTERPOLATED_INVALID = -1,
107a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_NONE    = 0,
108a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_PARTIAL,
109a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_MIXED,
110a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_FULL
111a7748839SVaclav Hapla } DMPlexInterpolatedFlag;
112a1cc8097SBarry Smith 
113a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
114a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *);
115a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolatePointSF(DM, PetscSF);
116a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolated(DM, DMPlexInterpolatedFlag *);
117a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolatedCollective(DM, DMPlexInterpolatedFlag *);
118a7748839SVaclav Hapla 
1198b51c812SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *);
120552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
121552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
122ef48cebcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *);
12308a22f4bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *);
12418e14f0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateLabelField(DM, DMLabel, Vec *);
125552f7358SJed Brown 
12675d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
1273ded2ed9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *);
12875d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
12975d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
130ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointDepth(DM, PetscInt, PetscInt *);
1310c0a32dcSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetPointHeight(DM, PetscInt, PetscInt *);
132ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellTypeLabel(DM, DMLabel *);
133ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellType(DM, PetscInt, DMPolytopeType *);
134412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetCellType(DM, PetscInt, DMPolytopeType);
135ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellTypes(DM);
13696ca5757SLisandro Dalcin PETSC_EXTERN PetscErrorCode DMPlexInvertCell(DMPolytopeType, PetscInt[]);
13796ca5757SLisandro Dalcin PETSC_EXTERN PetscErrorCode DMPlexReorderCell(DM, PetscInt, PetscInt[]);
13896ca5757SLisandro Dalcin 
13975d3a19aSMatthew G. Knepley /* Topological Operations */
140552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
141552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
142552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
143552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
144552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
145552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
146552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
147552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
14807218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
149b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
150552f7358SJed Brown 
151b7f5c055SJed Brown /*E
15287497f52SBarry Smith    DMPlexTPSType - Type of triply-periodic surface for a `DMPLEX`
153b7f5c055SJed Brown 
15416a05f60SBarry Smith    Values:
15516a05f60SBarry Smith +  `DMPLEX_TPS_SCHWARZ_P` - Schwarz Primitive surface, defined by the equation cos(x) + cos(y) + cos(z) = 0.
15616a05f60SBarry Smith -  `DMPLEX_TPS_GYROID`    - Gyroid surface, defined by the equation sin(x)cos(y) + sin(y)cos(z) + sin(z)cos(x) = 0
157b7f5c055SJed Brown 
158b7f5c055SJed Brown    Level: intermediate
159b7f5c055SJed Brown 
160b7f5c055SJed Brown    Developer Note:
161b7f5c055SJed Brown    Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
162b7f5c055SJed Brown 
16316a05f60SBarry Smith .seealso: `DMPLEX`, `DMPlexCreateTPSMesh()`
164b7f5c055SJed Brown E*/
165b7f5c055SJed Brown typedef enum {
166b7f5c055SJed Brown   DMPLEX_TPS_SCHWARZ_P,
167b7f5c055SJed Brown   DMPLEX_TPS_GYROID
168b7f5c055SJed Brown } DMPlexTPSType;
169b7f5c055SJed Brown PETSC_EXTERN const char *const DMPlexTPSTypes[];
170b7f5c055SJed Brown 
17175d3a19aSMatthew G. Knepley /* Mesh Generation */
172552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char[], PetscBool, DM *);
173930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM);
1749371c9d4SSatish Balay 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[]));
1750fdc7489SMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscReal, DM *);
176768d5fceSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
1779318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscBool, DM *);
17851a74b61SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, PetscReal, DM *);
17951a74b61SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBallMesh(MPI_Comm, PetscInt, PetscReal, DM *);
1809318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, DMBoundaryType, DM *);
1811436d7faSJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateTPSMesh(MPI_Comm, DMPlexTPSType, const PetscInt[], const DMBoundaryType[], PetscBool, PetscInt, PetscInt, PetscReal, DM *);
18224119c2aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *);
18300dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *);
184*1fcf445aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, PetscBool, const PetscReal[], const PetscReal[], DM *);
185552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
186c1cad2e7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModel(DM);
187068a5610SStefano Zampini 
1885f06a3ddSJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM, PetscSF);
1895f06a3ddSJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM, PetscSF *);
1905f06a3ddSJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM, const PetscScalar[]);
1914e2e9504SJed Brown 
1927f9d8d6cSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheck(DM);
193ca8062c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM);
19425c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscInt);
19525c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscInt);
196bb6a34a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckGeometry(DM);
197d7d32a9aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM, PetscSF, PetscBool);
198a8432d5bSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckInterfaceCones(DM);
19943fa8764SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckCellShape(DM, PetscBool, PetscReal);
200f108dbd7SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexComputeOrthogonalQuality(DM, PetscFV, PetscReal, Vec *, DMLabel *);
20175d3a19aSMatthew G. Knepley 
202d9deefdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *);
203776fd405SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *);
204d9deefdfSMatthew G. Knepley 
205cd7e8a5eSksagiyam PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], const char[], PetscBool, DM *);
206d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
207d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char[], PetscBool, DM *);
208d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *);
209d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *);
210d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *);
211d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char[], PetscBool, DM *);
212d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *);
213d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char[], PetscBool, DM *);
214d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char[], PetscBool, DM *);
215f2801cd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char[], PetscBool, DM *);
2167bee2925SMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm, const char[], DM *);
217c1cad2e7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSLiteFromFile(MPI_Comm, const char[], DM *);
218d6218766SMatthew G. Knepley 
2196823f3c5SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo);
2201e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetId(PetscViewer, int *);
2216823f3c5SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer, PetscInt);
2226823f3c5SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer, PetscInt *);
2231e50132fSMatthew G. Knepley 
22475d3a19aSMatthew G. Knepley /* Mesh Partitioning and Distribution */
225d2522c19Sksagiyam #define DMPLEX_OVERLAP_MANUAL -1
226d2522c19Sksagiyam 
227552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
2285680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *);
22971bb2955SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner);
2303fa7752dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *);
2311fd9873aSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel);
2325abbe4feSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel);
23324d039d7SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel);
234be200f8dSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel);
235aa3148a8SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *);
23698ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool);
23798ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *);
2385fa78c88SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsDistributed(DM, PetscBool *);
23980cf41d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF *, DM *);
240b9f40539SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *);
241cb54e036SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetOverlap(DM, PetscInt *);
242c506a872SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetOverlap(DM, DM, PetscInt);
243e600fa54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeGetDefault(DM, PetscBool *);
244e600fa54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeSetDefault(DM, PetscBool);
245552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM, PetscSF, PetscSection, Vec, PetscSection, Vec);
246b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *);
24793d501b3SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM, PetscSF, PetscSection, MPI_Datatype, void *, PetscSection, void **) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(5, 4) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(7, 4);
2488b879b41SFlorian Wechsung PETSC_EXTERN PetscErrorCode DMPlexRebalanceSharedPoints(DM, PetscInt, PetscBool, PetscBool, PetscBool *);
24945800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM);
250a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetGatherDM(DM, PetscSF *, DM *);
251a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, PetscSF *, DM *);
2523c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM, PetscInt, PetscInt *, PetscInt[], void *), void *);
2533c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM, PetscInt, PetscInt *, PetscInt[], void *), void **);
254a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM, PetscBool);
255a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM, PetscBool *);
25670034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]);
257f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF);
258f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *);
2591d1f2f2aSksagiyam PETSC_EXTERN PetscErrorCode DMPlexDistributionSetName(DM, const char[]);
2601d1f2f2aSksagiyam PETSC_EXTERN PetscErrorCode DMPlexDistributionGetName(DM, const char *[]);
26175d3a19aSMatthew G. Knepley 
2626bc1bd01Sksagiyam /*E
26387497f52SBarry Smith    DMPlexReorderDefaultFlag - Flag indicating whether the `DMPLEX` should be reordered by default
2646bc1bd01Sksagiyam 
26516a05f60SBarry Smith    Values:
26616a05f60SBarry Smith +  `DMPLEX_REORDER_DEFAULT_NOTSET` - Flag not set.
26716a05f60SBarry Smith .  `DMPLEX_REORDER_DEFAULT_FALSE`  - Do not reorder by default.
26816a05f60SBarry Smith -  `DMPLEX_REORDER_DEFAULT_TRUE`   - Reorder by default.
2696bc1bd01Sksagiyam 
2706bc1bd01Sksagiyam    Level: intermediate
2716bc1bd01Sksagiyam 
2726bc1bd01Sksagiyam    Developer Note:
2736bc1bd01Sksagiyam    Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
2746bc1bd01Sksagiyam 
27587497f52SBarry Smith    Could be replaced with `PETSC_BOOL3`
27687497f52SBarry Smith 
27787497f52SBarry Smith .seealso: `DMPlexReorderSetDefault()`, `DMPlexReorderGetDefault()`, `DMPlexGetOrdering()`, `DMPlexPermute()`
2786bc1bd01Sksagiyam E*/
2796bc1bd01Sksagiyam typedef enum {
2806bc1bd01Sksagiyam   DMPLEX_REORDER_DEFAULT_NOTSET = -1,
2816bc1bd01Sksagiyam   DMPLEX_REORDER_DEFAULT_FALSE  = 0,
2826bc1bd01Sksagiyam   DMPLEX_REORDER_DEFAULT_TRUE
2836bc1bd01Sksagiyam } DMPlexReorderDefaultFlag;
284c99f2573SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *);
2856913077dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering1D(DM, IS *);
286f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *);
2876bc1bd01Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexReorderGetDefault(DM, DMPlexReorderDefaultFlag *);
2886bc1bd01Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexReorderSetDefault(DM, DMPlexReorderDefaultFlag);
289f5cedc29SMatthew G. Knepley 
290b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *);
291b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *);
292b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *);
293c506a872SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointSF(DM, PetscSF, PetscBool, PetscSF *);
294dccdeccaSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabel(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *);
295c506a872SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabelFromLabels(DM, PetscInt, const DMLabel[], const PetscInt[], PetscInt, const DMLabel[], const PetscInt[], PetscSection, IS, PetscSection, IS, DMLabel *);
29646f9b1c3SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *);
29745800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *);
298b0a623aaSMatthew G. Knepley 
29975d3a19aSMatthew G. Knepley /* Submesh Support */
300158acfadSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM *);
301caf9e14dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, PetscInt, DMLabel *, DMLabel *, DM *, DM *);
302552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel *);
303552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
30497d8846cSMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpointIS(DM, IS *);
305a6e0b375SMatthew G. Knepley 
306a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetEnclosureRelation(DM, DM, DMEnclosureType *);
307a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetEnclosurePoint(DM, DM, DMEnclosureType, PetscInt, PetscInt *);
308552f7358SJed Brown 
3092be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
310caf9e14dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscInt, PetscBool, DM);
3116cf0e42fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel);
312a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddFaceCells(DM, DMLabel);
313f402d5e4SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel);
314552f7358SJed Brown 
31575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
31675d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
31775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
31875d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
3191e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal[], PetscReal *));
3201e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal[], PetscReal *));
3212389894bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *);
3220aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *);
3230aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool);
324552f7358SJed Brown 
32518ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */
32618ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
327f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt[], PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscBool *);
32818ad9376SMatthew G. Knepley 
3290f0bf202SMatthew G. Knepley /* Geometry Support */
3300f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *);
3310f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal);
332741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]);
333741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]);
334741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]);
3350f0bf202SMatthew G. Knepley 
336c4eade1cSMatthew G. Knepley /* Point Location */
337c4eade1cSMatthew G. Knepley typedef struct _PetscGridHash *PetscGridHash;
338c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *);
339c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscGridHashEnlarge(PetscGridHash, const PetscScalar[]);
340c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscGridHashSetGrid(PetscGridHash, const PetscInt[], const PetscReal[]);
3411c6dfc3eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar[], PetscInt[], PetscInt[]);
342c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscGridHashDestroy(PetscGridHash *);
343d3e1f4ccSVaclav Hapla PETSC_EXTERN PetscErrorCode    DMPlexFindVertices(DM, Vec, PetscReal, IS *);
344c4eade1cSMatthew G. Knepley 
345834e62ceSMatthew G. Knepley /* FVM Support */
346cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal[], PetscReal[]);
3470f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *);
348856ac710SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *);
349b27d5b9eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *);
350834e62ceSMatthew G. Knepley 
351552f7358SJed Brown /* FEM Support */
3523e9753d6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
3533e9753d6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGradientDM(DM, PetscFV, DM *);
354bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
35556cf3b9cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
3569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec);
3579371c9d4SSatish Balay 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);
3589371c9d4SSatish Balay 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);
3599371c9d4SSatish Balay 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);
360e752be1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel);
361d7ddef95SMatthew G. Knepley 
362baef929fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, DMLabel[], const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const IS[], const IS[], IS, PetscSection *);
363be36d101SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection *);
36475d3a19aSMatthew G. Knepley 
3658e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
366dfccc68fSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
3676858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellCoordinates(DM, PetscInt, PetscBool *, PetscInt *, const PetscScalar *[], PetscScalar *[]);
3686858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCellCoordinates(DM, PetscInt, PetscBool *, PetscInt *, const PetscScalar *[], PetscScalar *[]);
369d6143a4eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[]);
3709d150b73SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[]);
3713ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexShearGeometry(DM, DMDirection, PetscReal[]);
3729371c9d4SSatish Balay 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[]));
373d6143a4eSToby Isaac 
3747c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
3757c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
376552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
3777c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
37871f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureGeneral(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
37971f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
38071f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
3810405ed22SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
3827c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]);
383ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection);
384bc1eb3faSJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetClosurePermutationTensor(DM, PetscInt, PetscSection);
385552f7358SJed Brown 
386552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char[], PetscInt *, DM *);
3877db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *);
3880e18dc48SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexReorderCohesiveSupports(DM);
389552f7358SJed Brown 
390552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
391552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
392552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
393e6139122SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGhostCellStratum(DM, PetscInt *, PetscInt *);
394412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSimplexOrBoxCells(DM, PetscInt, PetscInt *, PetscInt *);
3959318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexIsSimplex(DM, PetscBool *);
396552f7358SJed Brown 
3972f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
3982f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
3992f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
4002f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
4012f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
4022f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
4032f856554SMatthew G. Knepley 
404552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
405552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
406552f7358SJed Brown 
407552f7358SJed Brown typedef struct {
408552f7358SJed Brown   DM    dm;
40928d58a37SPierre Jolivet   Vec   u; /* The base vector for the Jacobian action J(u) x */
410552f7358SJed Brown   Mat   J; /* Preconditioner for testing */
411552f7358SJed Brown   void *user;
412552f7358SJed Brown } JacActionCtx;
413552f7358SJed Brown 
414b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt);
415b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt *);
4161b32699bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetActivePoint(DM, PetscInt *);
4171b32699bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetActivePoint(DM, PetscInt);
4189371c9d4SSatish Balay 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);
419574a98acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
4200163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]);
4210163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, Vec);
422338f77d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *);
423b8feb594SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *);
4249371c9d4SSatish Balay 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 *);
425cf51de39SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, PetscBool, Mat, void *);
42668132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *);
4275f0b18bfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeClementInterpolant(DM, Vec, Vec);
4281555c271SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec);
4297c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *);
4302cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *);
4312cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *);
4321c41a8caSMatthew G. Knepley 
43356cf3b9cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, PetscInt, MatNullSpace *);
434fd86e548SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *);
435026175e5SToby Isaac 
4369f520fc2SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM, void *, void *, void *);
437bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *);
4380f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *);
4390f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *);
4400c290148SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, PetscWeakForm, PetscFormKey, Vec, Vec, Vec);
44145480ffeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, PetscWeakForm, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat);
4420f2d7e86SMatthew G. Knepley 
443ef68eab9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *);
444924a1b8fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
445adbe6fbbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *);
446756a1f44SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
447d64986aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM, PetscReal, Vec, Vec, void *);
448adbe6fbbSMatthew G. Knepley 
4491c41a8caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
4505d16530eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM, Vec, Vec);
451a68b90caSToby Isaac 
452a17985deSToby Isaac /* anchors */
453a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection *, IS *);
454a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS);
455d6a7ad0dSToby Isaac /* tree */
456d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM);
457d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM *);
458dcbd3bf7SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
459da43764aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM *);
460b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt[], PetscInt[]);
461b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]);
462d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *);
463d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]);
4646f5f1567SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *);
4653b490a17SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *);
466bbcf1e96SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM, Vec, DM, Vec, PetscSF, PetscSF, PetscInt *, PetscInt *, PetscBool, PetscReal);
467fa534816SMatthew G. Knepley 
468c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMonitorThroughput(DM, void *);
469c0f0dcc3SMatthew G. Knepley 
470fa534816SMatthew G. Knepley /* natural order */
471fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *);
472f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF);
473f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *);
474fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec);
475fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec);
476fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec);
477fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec);
47809cf685aSAlexis Marboeuf PETSC_EXTERN PetscErrorCode DMPlexCreateNaturalVector(DM, Vec *);
4790bbe5d31Sbarral 
4800bbe5d31Sbarral /* mesh adaptation */
4810f7e110dSbarral PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char[], DM *);
482d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
4839bfa2b02SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetFromOptions(DM);
48431b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetIsotropic(DM, PetscBool);
48531b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIsIsotropic(DM, PetscBool *);
48693520041SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetUniform(DM, PetscBool);
48793520041SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIsUniform(DM, PetscBool *);
48831b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM, PetscBool);
48931b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM, PetscBool *);
490cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoInsertion(DM, PetscBool);
491cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoInsertion(DM, PetscBool *);
492cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSwapping(DM, PetscBool);
493cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoSwapping(DM, PetscBool *);
494cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoMovement(DM, PetscBool);
495cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoMovement(DM, PetscBool *);
49676f360caSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSurf(DM, PetscBool);
49776f360caSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoSurf(DM, PetscBool *);
49831b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM, PetscReal);
49931b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM, PetscReal *);
50031b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM, PetscReal);
50131b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM, PetscReal *);
50231b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM, PetscReal);
50331b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM, PetscReal *);
50431b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetTargetComplexity(DM, PetscReal);
50531b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetTargetComplexity(DM, PetscReal *);
50631b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNormalizationOrder(DM, PetscReal);
50731b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetNormalizationOrder(DM, PetscReal *);
508cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetGradationFactor(DM, PetscReal);
509cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetGradationFactor(DM, PetscReal *);
510ae8b063eSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetHausdorffNumber(DM, PetscReal);
511ae8b063eSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetHausdorffNumber(DM, PetscReal *);
512cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetVerbosity(DM, PetscInt);
513cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetVerbosity(DM, PetscInt *);
514cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNumIterations(DM, PetscInt);
515cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetNumIterations(DM, PetscInt *);
5160f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricCreate(DM, PetscInt, Vec *);
5170f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricCreateUniform(DM, PetscInt, PetscReal, Vec *);
5180f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricCreateIsotropic(DM, PetscInt, Vec, Vec *);
519f6db075bSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricDeterminantCreate(DM, PetscInt, Vec *, DM *);
5205241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricEnforceSPD(DM, Vec, PetscBool, PetscBool, Vec, Vec);
5215241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNormalize(DM, Vec, PetscBool, PetscBool, Vec, Vec);
5225241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricAverage(DM, PetscInt, PetscReal[], Vec[], Vec);
5235241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricAverage2(DM, Vec, Vec, Vec);
5245241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricAverage3(DM, Vec, Vec, Vec, Vec);
5255241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection(DM, PetscInt, Vec[], Vec);
5265241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection2(DM, Vec, Vec, Vec);
5275241a91fSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection3(DM, Vec, Vec, Vec, Vec);
528ca3d3a14SMatthew G. Knepley 
529ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToLocalBasis(DM, Vec);
530ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLocalToGlobalBasis(DM, Vec);
531ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBasisRotation(DM, PetscReal, PetscReal, PetscReal);
532412e9a14SMatthew G. Knepley 
533c9ad657eSksagiyam /* storage version */
53481da5c8cSVaclav Hapla #define DMPLEX_STORAGE_VERSION_FIRST  "1.0.0"
535c9ad657eSksagiyam #define DMPLEX_STORAGE_VERSION_STABLE "1.0.0"
5369ed1248dSVaclav Hapla #define DMPLEX_STORAGE_VERSION_LATEST "3.0.0"
537c9ad657eSksagiyam 
5387f96f51bSksagiyam PETSC_EXTERN PetscErrorCode DMPlexTopologyView(DM, PetscViewer);
53977b8e257Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexCoordinatesView(DM, PetscViewer);
540bd6565f1Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexLabelsView(DM, PetscViewer);
541021affd3Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexSectionView(DM, PetscViewer, DM);
5423e97647fSksagiyam PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorView(DM, PetscViewer, DM, Vec);
5433e97647fSksagiyam PETSC_EXTERN PetscErrorCode DMPlexLocalVectorView(DM, PetscViewer, DM, Vec);
5446913077dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexVecView1D(DM, PetscInt, Vec[], PetscViewer);
5457f96f51bSksagiyam 
546dec9e869Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexTopologyLoad(DM, PetscViewer, PetscSF *);
547c9ad657eSksagiyam PETSC_EXTERN PetscErrorCode DMPlexCoordinatesLoad(DM, PetscViewer, PetscSF);
548e6368b79SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexLabelsLoad(DM, PetscViewer, PetscSF);
549f84dd6b4Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexSectionLoad(DM, PetscViewer, DM, PetscSF, PetscSF *, PetscSF *);
5508be3dfe1Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorLoad(DM, PetscViewer, DM, PetscSF, Vec);
5518be3dfe1Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexLocalVectorLoad(DM, PetscViewer, DM, PetscSF, Vec);
552ea8e1828Sksagiyam 
553a2c9b50fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMPlexGetLocalOffsets(DM, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt **);
55452e7713aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetLocalOffsetsSupport(DM, DMLabel, PetscInt, PetscInt *, PetscInt *, PetscInt *, PetscInt **, PetscInt **);
555a2c9b50fSJeremy L Thompson 
556c50b2d26SMatthew G. Knepley /* point queue */
557c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueCreate(PetscInt, DMPlexPointQueue *);
558c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueDestroy(DMPlexPointQueue *);
559c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueEnsureSize(DMPlexPointQueue);
560c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueEnqueue(DMPlexPointQueue, PetscInt);
561c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueDequeue(DMPlexPointQueue, PetscInt *);
562c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueFront(DMPlexPointQueue, PetscInt *);
563c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueBack(DMPlexPointQueue, PetscInt *);
564c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscBool      DMPlexPointQueueEmpty(DMPlexPointQueue);
565c50b2d26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointQueueEmptyCollective(PetscObject, DMPlexPointQueue, PetscBool *);
566c50b2d26SMatthew G. Knepley 
567552f7358SJed Brown #endif
568