xref: /petsc/include/petscdmplex.h (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
1552f7358SJed Brown /*
2552f7358SJed Brown   DMPlex, for parallel unstructured distributed mesh problems.
3552f7358SJed Brown */
426bd1501SBarry Smith #if !defined(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 
183c41b853SStefano Zampini PETSC_EXTERN PetscErrorCode PetscPartitionerDMPlexPartition(PetscPartitioner, DM, PetscSection, PetscSection, IS *);
1977623264SMatthew G. Knepley 
205e488331SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellList(DM, PetscInt, PetscInt, PetscInt, const PetscInt[]);
21be8c289dSNicolas Barral PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellListParallel(DM, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscSF *, PetscInt **);
221edcf0b2SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM, PetscInt, const PetscReal[]);
231edcf0b2SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM, PetscInt, PetscSF, const PetscReal[]);
24b09969d6SVaclav Hapla 
25552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*);
2627c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *);
27a4a685f2SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], DM*);
28a4a685f2SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMPlexCreateFromCellListPetsc() (since version 3.14)") PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*);
29be8c289dSNicolas Barral PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], PetscSF *, PetscInt **, DM *);
30a4a685f2SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMPlexCreateFromCellListParallelPetsc() (since version 3.14)") PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const PetscReal[], PetscSF *, 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);
38feb68dd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt);
39552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
400ce7577fSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeTuple(DM, IS, PetscSection *, IS *);
41af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
42af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexRestoreConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
43af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursiveVertices(DM, IS, IS *);
44552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
459f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
4677c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
47552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
48552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
49552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
50552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
51552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
52552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
53552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
54552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
558cb4d582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *);
56552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
57552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
58552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
59552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
614e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*);
62b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrientPoint(DM,PetscInt,PetscInt);
63d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM);
64e101f074SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
65552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*);
66081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,void*);
67a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*);
68a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*);
69a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*);
70081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,void*);
71552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*);
72552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*);
73a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*);
74a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*);
75a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*);
76081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*, void*);
77552f7358SJed Brown 
78a7748839SVaclav Hapla /* Topological interpolation */
79a7748839SVaclav Hapla PETSC_EXTERN const char * const DMPlexInterpolatedFlags[];
809ade3219SVaclav Hapla /*E
819ade3219SVaclav Hapla    DMPlexInterpolatedFlag - Describes level of topological interpolatedness.
829ade3219SVaclav Hapla      It is a local or collective property depending on whether it is returned by DMPlexIsInterpolated() or DMPlexIsInterpolatedCollective().
839ade3219SVaclav Hapla 
849ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_INVALID - Uninitialized value (internal use only; never returned by DMPlexIsInterpolated() or DMPlexIsInterpolatedCollective())
859ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_NONE    - Mesh is not interpolated
869ade3219SVaclav Hapla $  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)
879ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_MIXED   - Can be returned only by DMPlexIsInterpolatedCollective(), meaning that DMPlexIsInterpolated() returns different interpolatedness on different ranks
889ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_FULL    - Mesh is fully interpolated
899ade3219SVaclav Hapla 
909ade3219SVaclav Hapla    Level: intermediate
919ade3219SVaclav Hapla 
929ade3219SVaclav Hapla    Developer Note:
939ade3219SVaclav Hapla    Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
949ade3219SVaclav Hapla 
959ade3219SVaclav Hapla .seealso: DMPlexIsInterpolated(), DMPlexIsInterpolatedCollective(), DMPlexInterpolate(), DMPlexUninterpolate()
969ade3219SVaclav Hapla E*/
97a7748839SVaclav Hapla typedef enum {
987d0f5628SVaclav Hapla   DMPLEX_INTERPOLATED_INVALID = -1,
99a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_NONE = 0,
100a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_PARTIAL,
101a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_MIXED,
102a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_FULL
103a7748839SVaclav Hapla } DMPlexInterpolatedFlag;
104a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
105a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *);
106a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolatePointSF(DM, PetscSF);
107a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolated(DM, DMPlexInterpolatedFlag *);
108a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolatedCollective(DM, DMPlexInterpolatedFlag *);
109a7748839SVaclav Hapla 
1108b51c812SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *);
111552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
112552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
113ef48cebcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *);
11408a22f4bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *);
11518e14f0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateLabelField(DM, DMLabel, Vec *);
116552f7358SJed Brown 
11775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
1183ded2ed9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *);
11975d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
12075d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
121ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointDepth(DM, PetscInt, PetscInt *);
1220c0a32dcSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetPointHeight(DM, PetscInt, PetscInt *);
123ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellTypeLabel(DM, DMLabel *);
124ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellType(DM, PetscInt, DMPolytopeType *);
125412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetCellType(DM, PetscInt, DMPolytopeType);
126ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellTypes(DM);
12796ca5757SLisandro Dalcin PETSC_EXTERN PetscErrorCode DMPlexInvertCell(DMPolytopeType, PetscInt[]);
12896ca5757SLisandro Dalcin PETSC_EXTERN PetscErrorCode DMPlexReorderCell(DM, PetscInt, PetscInt[]);
12996ca5757SLisandro Dalcin 
13075d3a19aSMatthew G. Knepley /* Topological Operations */
131552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
132552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
133552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
134552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
135552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
136552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
137552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
138552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
139b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
140b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
141552f7358SJed Brown 
14275d3a19aSMatthew G. Knepley /* Mesh Generation */
143552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *);
144930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM);
14548d16a33SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoordinateSpace(DM, PetscInt, void (*)(PetscInt, PetscInt, PetscInt,
14648d16a33SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
14748d16a33SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
14848d16a33SMatthew G. Knepley                                                         PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
1490fdc7489SMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscReal, DM *);
150768d5fceSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
1519318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscBool, DM *);
15251a74b61SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, PetscReal, DM *);
15351a74b61SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBallMesh(MPI_Comm, PetscInt, PetscReal, DM *);
1549318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, DMBoundaryType, DM *);
15524119c2aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *);
15600dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *);
157d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, const PetscReal[], const PetscReal[], DM *);
158552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
159c1cad2e7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModel(DM);
160068a5610SStefano Zampini 
161ca8062c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM);
16225c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscInt);
16325c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscInt);
164bb6a34a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckGeometry(DM);
16503da9461SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM);
166a8432d5bSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckInterfaceCones(DM);
16743fa8764SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckCellShape(DM, PetscBool, PetscReal);
168f108dbd7SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexComputeOrthogonalQuality(DM, PetscFV, PetscReal, Vec *, DMLabel *);
16975d3a19aSMatthew G. Knepley 
170d9deefdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *);
171776fd405SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *);
172d9deefdfSMatthew G. Knepley 
173cd7e8a5eSksagiyam PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], const char[], PetscBool, DM *);
174d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
175d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char [], PetscBool, DM *);
176d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *);
177d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *);
178d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *);
179d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char [], PetscBool, DM *);
180d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *);
181d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char [], PetscBool, DM *);
182d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char [], PetscBool, DM *);
183f2801cd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char [], PetscBool, DM *);
1847bee2925SMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm, const char [], DM *);
185c1cad2e7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSLiteFromFile(MPI_Comm, const char [], DM *);
186d6218766SMatthew G. Knepley 
1876823f3c5SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo);
1881e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetId(PetscViewer, int *);
1896823f3c5SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer,PetscInt);
1906823f3c5SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer,PetscInt*);
1911e50132fSMatthew G. Knepley 
19275d3a19aSMatthew G. Knepley /* Mesh Partitioning and Distribution */
193552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
1945680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *);
19571bb2955SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner);
1963fa7752dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *);
1971fd9873aSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel);
1985abbe4feSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel);
19924d039d7SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel);
200be200f8dSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel);
201aa3148a8SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *);
20298ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool);
20398ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *);
2045fa78c88SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsDistributed(DM, PetscBool *);
20580cf41d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF*, DM*);
206b9f40539SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *);
207cb54e036SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetOverlap(DM, PetscInt *);
208*e600fa54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeGetDefault(DM, PetscBool *);
209*e600fa54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeSetDefault(DM, PetscBool);
210552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM,PetscSF,PetscSection,Vec,PetscSection,Vec);
211b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *);
212552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM,PetscSF,PetscSection,MPI_Datatype,void*,PetscSection,void**);
2138b879b41SFlorian Wechsung PETSC_EXTERN PetscErrorCode DMPlexRebalanceSharedPoints(DM, PetscInt, PetscBool, PetscBool, PetscBool*);
21445800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM);
215a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetGatherDM(DM, PetscSF*, DM*);
216a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, PetscSF*, DM*);
2173c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM,PetscInt,PetscInt*,PetscInt[],void*),void*);
2183c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM,PetscInt,PetscInt*,PetscInt[],void*),void**);
219a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM,PetscBool);
220a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM,PetscBool*);
22170034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]);
222f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF);
223f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *);
22475d3a19aSMatthew G. Knepley 
225c99f2573SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *);
226f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *);
227f5cedc29SMatthew G. Knepley 
228b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *);
229b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *);
230b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *);
231dccdeccaSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabel(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *);
23246f9b1c3SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *);
23345800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *);
234b0a623aaSMatthew G. Knepley 
23575d3a19aSMatthew G. Knepley /* Submesh Support */
236158acfadSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM*);
2377db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, DMLabel *, DMLabel *, DM *, DM *);
238552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*);
239552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
24097d8846cSMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpointIS(DM, IS *);
241a6e0b375SMatthew G. Knepley 
242a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetEnclosureRelation(DM, DM, DMEnclosureType *);
243a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetEnclosurePoint(DM, DM, DMEnclosureType, PetscInt, PetscInt *);
244552f7358SJed Brown 
2452be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
2460f66a230SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscBool, DM);
2476cf0e42fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel);
248a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddFaceCells(DM, DMLabel);
249f402d5e4SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel);
250552f7358SJed Brown 
25175d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
25275d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
25375d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
25475d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
2551e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal [], PetscReal *));
2561e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal [], PetscReal *));
2572389894bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *);
2580aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *);
2590aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool);
260552f7358SJed Brown 
26118ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */
26218ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
263f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *);
26418ad9376SMatthew G. Knepley 
2650f0bf202SMatthew G. Knepley /* Geometry Support */
2660f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *);
2670f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal);
268741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]);
269741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]);
270741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]);
2710f0bf202SMatthew G. Knepley 
272c4eade1cSMatthew G. Knepley /* Point Location */
273c4eade1cSMatthew G. Knepley typedef struct _PetscGridHash *PetscGridHash;
274c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *);
275c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar []);
276c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt [], const PetscReal []);
2771c6dfc3eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar [], PetscInt [], PetscInt []);
278c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *);
279d3e1f4ccSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexFindVertices(DM, Vec, PetscReal, IS *);
280c4eade1cSMatthew G. Knepley 
281834e62ceSMatthew G. Knepley /* FVM Support */
282cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []);
2830f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *);
284856ac710SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *);
285b27d5b9eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *);
286834e62ceSMatthew G. Knepley 
287552f7358SJed Brown /* FEM Support */
288cf0b7c11SKarl Rupp PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *);
2893e9753d6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
2903e9753d6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGradientDM(DM,PetscFV,DM*);
291bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
29256cf3b9cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
2931c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[],
294b278463cSMatthew G. Knepley                                                                 PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec);
2951c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
296b278463cSMatthew G. Knepley                                                                      void (*)(PetscInt, PetscInt, PetscInt,
297b278463cSMatthew G. Knepley                                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
298b278463cSMatthew G. Knepley                                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
299191494d9SMatthew G. Knepley                                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[],
300b278463cSMatthew G. Knepley                                                                               PetscScalar[]),
301b278463cSMatthew G. Knepley                                                                      void *, Vec);
302ece3a9fcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
303ece3a9fcSMatthew G. Knepley                                                                        void (*)(PetscInt, PetscInt, PetscInt,
304ece3a9fcSMatthew G. Knepley                                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
305ece3a9fcSMatthew G. Knepley                                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
306ece3a9fcSMatthew G. Knepley                                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[],
307ece3a9fcSMatthew G. Knepley                                                                                 PetscScalar[]),
308ece3a9fcSMatthew G. Knepley                                                                        void *, Vec);
3091c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM, PetscReal, Vec, Vec, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
310b278463cSMatthew G. Knepley                                                               PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *, Vec);
311e752be1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel);
312d7ddef95SMatthew G. Knepley 
313baef929fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, DMLabel [], const PetscInt [], const PetscInt [], PetscInt, const PetscInt [], const IS [], const IS [], IS, PetscSection *);
314be36d101SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection*);
31575d3a19aSMatthew G. Knepley 
3168e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
317dfccc68fSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
318d6143a4eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]);
3199d150b73SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]);
3203ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexShearGeometry(DM, DMDirection, PetscReal[]);
3210139fca9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRemapGeometry(DM, PetscReal,
3220139fca9SMatthew G. Knepley                                                 void (*)(PetscInt, PetscInt, PetscInt,
3230139fca9SMatthew G. Knepley                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3240139fca9SMatthew G. Knepley                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3250139fca9SMatthew G. Knepley                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
326d6143a4eSToby Isaac 
3277c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
3287c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
329552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
3307c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
33171f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureGeneral(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
33271f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
33371f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
3340405ed22SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
3357c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]);
336ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection);
337bc1eb3faSJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetClosurePermutationTensor(DM, PetscInt, PetscSection);
338552f7358SJed Brown 
339552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *);
3407db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *);
341552f7358SJed Brown 
342552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
343552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
344552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
345e6139122SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGhostCellStratum(DM, PetscInt *, PetscInt *);
346412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSimplexOrBoxCells(DM, PetscInt, PetscInt *, PetscInt *);
3479318fe57SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexIsSimplex(DM, PetscBool *);
348552f7358SJed Brown 
3492f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
3502f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
3512f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
3522f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
3532f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
3542f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
3552f856554SMatthew G. Knepley 
356552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
357552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
358552f7358SJed Brown 
359552f7358SJed Brown typedef struct {
360552f7358SJed Brown   DM    dm;
36128d58a37SPierre Jolivet   Vec   u; /* The base vector for the Jacobian action J(u) x */
362552f7358SJed Brown   Mat   J; /* Preconditioner for testing */
363552f7358SJed Brown   void *user;
364552f7358SJed Brown } JacActionCtx;
365552f7358SJed Brown 
366b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt);
367b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt*);
3681b32699bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetActivePoint(DM, PetscInt*);
3691b32699bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetActivePoint(DM, PetscInt);
370e29c0834SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec,
371e29c0834SMatthew G. Knepley                                                     void (**)(PetscInt, PetscInt, PetscInt,
372e29c0834SMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
373e29c0834SMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
374e29c0834SMatthew G. Knepley                                                               PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec);
375574a98acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
3760163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]);
3770163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, Vec);
378338f77d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *);
379b8feb594SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *);
3809b6f715bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdIntegral(DM, Vec, DMLabel, PetscInt, const PetscInt[],
3819b6f715bSMatthew G. Knepley                                                     void (*)(PetscInt, PetscInt, PetscInt,
3829b6f715bSMatthew G. Knepley                                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3839b6f715bSMatthew G. Knepley                                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3849b6f715bSMatthew G. Knepley                                                              PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
3859b6f715bSMatthew G. Knepley                                                     PetscScalar *, void *);
386cf51de39SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, PetscBool, Mat, void *);
38768132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *);
3881555c271SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec);
3897c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *);
3902cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *);
3912cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *);
3921c41a8caSMatthew G. Knepley 
39356cf3b9cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, PetscInt, MatNullSpace *);
394fd86e548SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *);
395026175e5SToby Isaac 
3969f520fc2SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM,void *,void *,void *);
397bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *);
3980f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *);
3990f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *);
4000c290148SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, PetscWeakForm, PetscFormKey, Vec, Vec, Vec);
40145480ffeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, PetscWeakForm, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat);
4020f2d7e86SMatthew G. Knepley 
403ef68eab9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *);
404924a1b8fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
405adbe6fbbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *);
406756a1f44SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
407d64986aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM, PetscReal, Vec, Vec, void *);
408adbe6fbbSMatthew G. Knepley 
4091c41a8caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
4105d16530eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM,Vec,Vec);
411a68b90caSToby Isaac 
412a17985deSToby Isaac /* anchors */
413a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection*, IS*);
414a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS);
415d6a7ad0dSToby Isaac /* tree */
416d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM);
417d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM*);
418dcbd3bf7SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
419da43764aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM*);
420b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt [], PetscInt []);
421b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]);
422d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *);
423d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]);
4246f5f1567SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *);
4253b490a17SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *);
426bbcf1e96SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM,Vec,DM,Vec,PetscSF,PetscSF,PetscInt *,PetscInt *,PetscBool,PetscReal);
427fa534816SMatthew G. Knepley 
428c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMonitorThroughput(DM, void *);
429c0f0dcc3SMatthew G. Knepley 
430fa534816SMatthew G. Knepley /* natural order */
431fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *);
432f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF);
433f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *);
434fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec);
435fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec);
436fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec);
437fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec);
4380bbe5d31Sbarral 
4390bbe5d31Sbarral /* mesh adaptation */
4400f7e110dSbarral PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char [], DM *);
441d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
4429bfa2b02SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetFromOptions(DM);
44331b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetIsotropic(DM, PetscBool);
44431b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIsIsotropic(DM, PetscBool *);
44593520041SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetUniform(DM, PetscBool);
44693520041SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIsUniform(DM, PetscBool *);
44731b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM, PetscBool);
44831b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM, PetscBool *);
449cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoInsertion(DM, PetscBool);
450cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoInsertion(DM, PetscBool *);
451cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSwapping(DM, PetscBool);
452cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoSwapping(DM, PetscBool *);
453cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoMovement(DM, PetscBool);
454cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNoMovement(DM, PetscBool *);
45531b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM, PetscReal);
45631b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM, PetscReal *);
45731b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM, PetscReal);
45831b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM, PetscReal *);
45931b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM, PetscReal);
46031b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM, PetscReal *);
46131b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetTargetComplexity(DM, PetscReal);
46231b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetTargetComplexity(DM, PetscReal *);
46331b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNormalizationOrder(DM, PetscReal);
46431b70465SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetNormalizationOrder(DM, PetscReal *);
465cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetGradationFactor(DM, PetscReal);
466cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetGradationFactor(DM, PetscReal *);
467cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetVerbosity(DM, PetscInt);
468cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetVerbosity(DM, PetscInt *);
469cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricSetNumIterations(DM, PetscInt);
470cc2eee55SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricGetNumIterations(DM, PetscInt *);
4710f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricCreate(DM, PetscInt, Vec *);
4720f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricCreateUniform(DM, PetscInt, PetscReal, Vec *);
4730f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricCreateIsotropic(DM, PetscInt, Vec, Vec *);
4743f07679eSJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricEnforceSPD(DM, Vec, PetscBool, PetscBool, Vec *, Vec *);
47516522936SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricNormalize(DM, Vec, PetscBool, PetscBool, Vec *);
4760f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricAverage(DM, PetscInt, PetscReal[], Vec[], Vec *);
4770f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricAverage2(DM, Vec, Vec, Vec *);
4780f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricAverage3(DM, Vec, Vec, Vec, Vec *);
4790f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection(DM, PetscInt, Vec[], Vec *);
4800f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection2(DM, Vec, Vec, Vec *);
4810f923b09SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection3(DM, Vec, Vec, Vec, Vec *);
482ca3d3a14SMatthew G. Knepley 
483ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToLocalBasis(DM, Vec);
484ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLocalToGlobalBasis(DM, Vec);
485ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBasisRotation(DM, PetscReal, PetscReal, PetscReal);
486412e9a14SMatthew G. Knepley 
487c9ad657eSksagiyam /* storage version */
488c9ad657eSksagiyam #define DMPLEX_STORAGE_VERSION_STABLE "1.0.0"
489c9ad657eSksagiyam #define DMPLEX_STORAGE_VERSION_LATEST "2.0.0"
490c9ad657eSksagiyam 
4917f96f51bSksagiyam PETSC_EXTERN PetscErrorCode DMPlexTopologyView(DM, PetscViewer);
49277b8e257Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexCoordinatesView(DM, PetscViewer);
493bd6565f1Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexLabelsView(DM, PetscViewer);
494021affd3Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexSectionView(DM,PetscViewer,DM);
4953e97647fSksagiyam PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorView(DM,PetscViewer,DM,Vec);
4963e97647fSksagiyam PETSC_EXTERN PetscErrorCode DMPlexLocalVectorView(DM,PetscViewer,DM,Vec);
4977f96f51bSksagiyam 
498dec9e869Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexTopologyLoad(DM, PetscViewer, PetscSF*);
499c9ad657eSksagiyam PETSC_EXTERN PetscErrorCode DMPlexCoordinatesLoad(DM, PetscViewer, PetscSF);
500e6368b79SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexLabelsLoad(DM, PetscViewer, PetscSF);
501f84dd6b4Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexSectionLoad(DM,PetscViewer,DM,PetscSF,PetscSF*,PetscSF*);
5028be3dfe1Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorLoad(DM,PetscViewer,DM,PetscSF,Vec);
5038be3dfe1Sksagiyam PETSC_EXTERN PetscErrorCode DMPlexLocalVectorLoad(DM,PetscViewer,DM,PetscSF,Vec);
504ea8e1828Sksagiyam 
505a2c9b50fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMPlexGetLocalOffsets(DM, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt **);
506a2c9b50fSJeremy L Thompson 
507552f7358SJed Brown #endif
508