xref: /petsc/include/petscdmplex.h (revision 3e9753d60b406ac2fec1249fa8fb9cfa31a5ff24)
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>
16552f7358SJed Brown 
173c41b853SStefano Zampini PETSC_EXTERN PetscErrorCode PetscPartitionerDMPlexPartition(PetscPartitioner, DM, PetscSection, PetscSection, IS *);
1877623264SMatthew G. Knepley 
19552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*);
2027c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *);
210fde5641SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*);
22fba955ccSBarry Smith PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const PetscReal[], PetscSF *, DM *);
2306bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []);
248415267dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, PetscInt, PetscBool, DM*);
25412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCellByType(MPI_Comm, DMPolytopeType, DM *);
26a9074c1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetOptionsPrefix(DM, const char []);
27552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *);
28552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt);
29552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *);
30552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt);
31feb68dd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt);
32552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
330ce7577fSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeTuple(DM, IS, PetscSection *, IS *);
34af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
35af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexRestoreConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
36af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursiveVertices(DM, IS, IS *);
37552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
389f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
3977c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
40552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
41552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
42552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
43552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
44552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
45552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
46552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
47552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
488cb4d582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *);
49552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
50552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
51552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
52552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
53552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
544e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*);
55cb11e54dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexReverseCell(DM, PetscInt);
56fa01033eSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexOrientCell(DM,PetscInt,PetscInt,const PetscInt[]);
5727d0e99aSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCompareOrientations(DM, PetscInt, PetscInt, const PetscInt [], PetscInt *, PetscBool *);
58d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM);
59e101f074SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*);
61081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,void*);
62a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*);
63a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*);
64a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*);
65081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,void*);
66552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*);
67552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*);
68a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*);
69a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*);
70a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*);
71081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*, void*);
72552f7358SJed Brown 
73a7748839SVaclav Hapla /* Topological interpolation */
74a7748839SVaclav Hapla PETSC_EXTERN const char * const DMPlexInterpolatedFlags[];
759ade3219SVaclav Hapla /*E
769ade3219SVaclav Hapla    DMPlexInterpolatedFlag - Describes level of topological interpolatedness.
779ade3219SVaclav Hapla      It is a local or collective property depending on whether it is returned by DMPlexIsInterpolated() or DMPlexIsInterpolatedCollective().
789ade3219SVaclav Hapla 
799ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_INVALID - Uninitialized value (internal use only; never returned by DMPlexIsInterpolated() or DMPlexIsInterpolatedCollective())
809ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_NONE    - Mesh is not interpolated
819ade3219SVaclav 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)
829ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_MIXED   - Can be returned only by DMPlexIsInterpolatedCollective(), meaning that DMPlexIsInterpolated() returns different interpolatedness on different ranks
839ade3219SVaclav Hapla $  DMPLEX_INTERPOLATED_FULL    - Mesh is fully interpolated
849ade3219SVaclav Hapla 
859ade3219SVaclav Hapla    Level: intermediate
869ade3219SVaclav Hapla 
879ade3219SVaclav Hapla    Developer Note:
889ade3219SVaclav Hapla    Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
899ade3219SVaclav Hapla 
909ade3219SVaclav Hapla .seealso: DMPlexIsInterpolated(), DMPlexIsInterpolatedCollective(), DMPlexInterpolate(), DMPlexUninterpolate()
919ade3219SVaclav Hapla E*/
92a7748839SVaclav Hapla typedef enum {
937d0f5628SVaclav Hapla   DMPLEX_INTERPOLATED_INVALID = -1,
94a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_NONE = 0,
95a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_PARTIAL,
96a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_MIXED,
97a7748839SVaclav Hapla   DMPLEX_INTERPOLATED_FULL
98a7748839SVaclav Hapla } DMPlexInterpolatedFlag;
99a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
100a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *);
101a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolatePointSF(DM, PetscSF);
102a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolated(DM, DMPlexInterpolatedFlag *);
103a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolatedCollective(DM, DMPlexInterpolatedFlag *);
104a7748839SVaclav Hapla 
1058b51c812SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *);
106552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
107552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
108ef48cebcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *);
10908a22f4bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *);
11018e14f0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateLabelField(DM, DMLabel, Vec *);
111552f7358SJed Brown 
11275d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
1133ded2ed9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *);
11475d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
11575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
116ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointDepth(DM, PetscInt, PetscInt *);
1170c0a32dcSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetPointHeight(DM, PetscInt, PetscInt *);
118ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellTypeLabel(DM, DMLabel *);
119ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellType(DM, PetscInt, DMPolytopeType *);
120412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetCellType(DM, PetscInt, DMPolytopeType);
121ba2698f1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellTypes(DM);
12296ca5757SLisandro Dalcin PETSC_EXTERN PetscErrorCode DMPlexInvertCell(DMPolytopeType, PetscInt[]);
12396ca5757SLisandro Dalcin PETSC_EXTERN PetscErrorCode DMPlexReorderCell(DM, PetscInt, PetscInt[]);
12496ca5757SLisandro Dalcin 
12575d3a19aSMatthew G. Knepley 
12675d3a19aSMatthew G. Knepley /* Topological Operations */
127552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
128552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
129552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
130552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
131552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
132552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
133552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
134552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
135552f7358SJed Brown 
13675d3a19aSMatthew G. Knepley /* Mesh Generation */
137552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *);
1383a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegister(const char[],PetscErrorCode (*)(DM,PetscBool,DM*),PetscErrorCode (*)(DM,double*,DM*),PetscInt);
1393a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegisterAll(void);
1403a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegisterDestroy(void);
141930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM);
1421df5d5c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscBool, PetscReal, DM *);
14326492d91SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSquareBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []);
144552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []);
145768d5fceSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
14665a81367SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, DM *);
147dbc1dc17SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, PetscInt, DMBoundaryType, DM *);
14824119c2aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *);
14900dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *);
15000dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, DM *);
151552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
152068a5610SStefano Zampini 
153ca8062c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM);
15425c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscInt);
15525c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscInt);
156bb6a34a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckGeometry(DM);
15703da9461SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM);
158a8432d5bSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckInterfaceCones(DM);
15943fa8764SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckCellShape(DM, PetscBool, PetscReal);
160f108dbd7SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMPlexComputeOrthogonalQuality(DM, PetscFV, PetscReal, Vec *, DMLabel *);
16175d3a19aSMatthew G. Knepley 
162d9deefdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *);
163776fd405SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *);
164d9deefdfSMatthew G. Knepley 
165d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], PetscBool, DM *);
166d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
167d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char [], PetscBool, DM *);
168d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *);
169d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *);
170d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *);
171d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char [], PetscBool, DM *);
172d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *);
173d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char [], PetscBool, DM *);
174d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char [], PetscBool, DM *);
175f2801cd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char [], PetscBool, DM *);
176d6218766SMatthew G. Knepley 
1771e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetId(PetscViewer, int *);
1781e50132fSMatthew G. Knepley 
17975d3a19aSMatthew G. Knepley /* Mesh Partitioning and Distribution */
180552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
1815680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *);
18271bb2955SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner);
183ac17c9edSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *);
1843fa7752dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *);
185552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *);
1861fd9873aSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel);
1875abbe4feSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel);
18824d039d7SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel);
189be200f8dSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel);
190aa3148a8SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *);
19198ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool);
19298ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *);
1935fa78c88SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsDistributed(DM, PetscBool *);
19480cf41d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF*, DM*);
195b9f40539SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *);
196cb54e036SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetOverlap(DM, PetscInt *);
197552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM,PetscSF,PetscSection,Vec,PetscSection,Vec);
198b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *);
199552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM,PetscSF,PetscSection,MPI_Datatype,void*,PetscSection,void**);
2008b879b41SFlorian Wechsung PETSC_EXTERN PetscErrorCode DMPlexRebalanceSharedPoints(DM, PetscInt, PetscBool, PetscBool, PetscBool*);
20145800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM);
202a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetGatherDM(DM, PetscSF*, DM*);
203a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, PetscSF*, DM*);
2043c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM,PetscInt,PetscInt*,PetscInt[],void*),void*);
2053c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM,PetscInt,PetscInt*,PetscInt[],void*),void**);
20670034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseCone(DM,PetscBool);
20770034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseCone(DM,PetscBool*);
20870034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseClosure(DM,PetscBool);
20970034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseClosure(DM,PetscBool*);
210a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM,PetscBool);
211a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM,PetscBool*);
21270034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]);
213f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF);
214f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *);
21575d3a19aSMatthew G. Knepley 
216c99f2573SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *);
217f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *);
218f5cedc29SMatthew G. Knepley 
219b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *);
220b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *);
221b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *);
222dccdeccaSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabel(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *);
22346f9b1c3SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *);
22445800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *);
225b0a623aaSMatthew G. Knepley 
22675d3a19aSMatthew G. Knepley /* Submesh Support */
227158acfadSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM*);
2287db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, DMLabel *, DMLabel *, DM *, DM *);
229552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*);
230552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
23197d8846cSMatthew Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpointIS(DM, IS *);
232a6e0b375SMatthew G. Knepley 
233a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetEnclosureRelation(DM, DM, DMEnclosureType *);
234a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetEnclosurePoint(DM, DM, DMEnclosureType, PetscInt, PetscInt *);
235552f7358SJed Brown 
2362be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
2370f66a230SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscBool, DM);
2386cf0e42fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel);
239a6e0b375SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddFaceCells(DM, DMLabel);
240f402d5e4SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel);
241552f7358SJed Brown 
24275d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
24375d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
24475d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
24575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
2461e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal [], PetscReal *));
2471e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal [], PetscReal *));
2482389894bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *);
2490aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *);
2500aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool);
251412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellRefinerType(DM, DMPlexCellRefinerType *);
252412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetCellRefinerType(DM, DMPlexCellRefinerType);
253552f7358SJed Brown 
25418ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */
25518ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
256f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *);
25718ad9376SMatthew G. Knepley 
2580f0bf202SMatthew G. Knepley /* Geometry Support */
2590f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *);
2600f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal);
261741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]);
262741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]);
263741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]);
2640f0bf202SMatthew G. Knepley 
265c4eade1cSMatthew G. Knepley /* Point Location */
266c4eade1cSMatthew G. Knepley typedef struct _PetscGridHash *PetscGridHash;
267c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *);
268c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar []);
269c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt [], const PetscReal []);
2701c6dfc3eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar [], PetscInt [], PetscInt []);
271c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *);
2723985bb02SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexFindVertices(DM, PetscInt, const PetscReal [], PetscReal, PetscInt []);
273c4eade1cSMatthew G. Knepley 
274834e62ceSMatthew G. Knepley /* FVM Support */
275cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []);
2760f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *);
277856ac710SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *);
278b27d5b9eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *);
279834e62ceSMatthew G. Knepley 
280552f7358SJed Brown /* FEM Support */
281cf0b7c11SKarl Rupp PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *);
282*3e9753d6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGeometryFVM(DM,Vec*,Vec*,PetscReal*);
283*3e9753d6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGradientDM(DM,PetscFV,DM*);
284bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
2851c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[],
286b278463cSMatthew G. Knepley                                                                 PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec );
2871c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
288b278463cSMatthew G. Knepley                                                                      void (*)(PetscInt, PetscInt, PetscInt,
289b278463cSMatthew G. Knepley                                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
290b278463cSMatthew G. Knepley                                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
291191494d9SMatthew G. Knepley                                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[],
292b278463cSMatthew G. Knepley                                                                               PetscScalar[]),
293b278463cSMatthew G. Knepley                                                                      void *, Vec);
294ece3a9fcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
295ece3a9fcSMatthew G. Knepley                                                                        void (*)(PetscInt, PetscInt, PetscInt,
296ece3a9fcSMatthew G. Knepley                                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
297ece3a9fcSMatthew G. Knepley                                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
298ece3a9fcSMatthew G. Knepley                                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[],
299ece3a9fcSMatthew G. Knepley                                                                                 PetscScalar[]),
300ece3a9fcSMatthew G. Knepley                                                                        void *, Vec);
3011c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM, PetscReal, Vec, Vec, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
302b278463cSMatthew G. Knepley                                                               PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *, Vec);
303e752be1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel);
304d7ddef95SMatthew G. Knepley 
305baef929fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, DMLabel [], const PetscInt [], const PetscInt [], PetscInt, const PetscInt [], const IS [], const IS [], IS, PetscSection *);
306be36d101SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection*);
30775d3a19aSMatthew G. Knepley 
3088e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
309dfccc68fSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
310d6143a4eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]);
3119d150b73SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]);
3123ee9839eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexShearGeometry(DM, DMDirection, PetscReal[]);
3130139fca9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRemapGeometry(DM, PetscReal,
3140139fca9SMatthew G. Knepley                                                 void (*)(PetscInt, PetscInt, PetscInt,
3150139fca9SMatthew G. Knepley                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3160139fca9SMatthew G. Knepley                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3170139fca9SMatthew G. Knepley                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
318d6143a4eSToby Isaac 
3197c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
3207c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
321552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
3227c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
32371f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureGeneral(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
32471f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
32571f0bbf9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
3260405ed22SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
3277c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]);
328ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection);
329bc1eb3faSJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetClosurePermutationTensor(DM, PetscInt, PetscSection);
330552f7358SJed Brown 
331552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *);
3327db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *);
333552f7358SJed Brown 
334552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
335552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
336552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
337e6139122SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetGhostCellStratum(DM, PetscInt *, PetscInt *);
338412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSimplexOrBoxCells(DM, PetscInt, PetscInt *, PetscInt *);
339552f7358SJed Brown 
3402f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
3412f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
3422f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
3432f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
3442f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
3452f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
3462f856554SMatthew G. Knepley 
347552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
348552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
349552f7358SJed Brown 
350552f7358SJed Brown typedef struct {
351552f7358SJed Brown   DM    dm;
35228d58a37SPierre Jolivet   Vec   u; /* The base vector for the Jacobian action J(u) x */
353552f7358SJed Brown   Mat   J; /* Preconditioner for testing */
354552f7358SJed Brown   void *user;
355552f7358SJed Brown } JacActionCtx;
356552f7358SJed Brown 
3573351dd3dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM, Vec);
358b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt);
359b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt*);
3601b32699bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetActivePoint(DM, PetscInt*);
3611b32699bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetActivePoint(DM, PetscInt);
362e29c0834SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec,
363e29c0834SMatthew G. Knepley                                                     void (**)(PetscInt, PetscInt, PetscInt,
364e29c0834SMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
365e29c0834SMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
366e29c0834SMatthew G. Knepley                                                               PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec);
367574a98acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
3680163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]);
3690163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, Vec);
370338f77d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *);
371b8feb594SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *);
3729b6f715bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdIntegral(DM, Vec, DMLabel, PetscInt, const PetscInt[],
3739b6f715bSMatthew G. Knepley                                                     void (*)(PetscInt, PetscInt, PetscInt,
3749b6f715bSMatthew G. Knepley                                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3759b6f715bSMatthew G. Knepley                                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3769b6f715bSMatthew G. Knepley                                                              PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
3779b6f715bSMatthew G. Knepley                                                     PetscScalar *, void *);
378cf51de39SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, PetscBool, Mat, void *);
37968132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *);
3801555c271SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec);
3817c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *);
3822cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *);
3832cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *);
3841c41a8caSMatthew G. Knepley 
385026175e5SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, MatNullSpace *);
386fd86e548SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *);
387026175e5SToby Isaac 
3889f520fc2SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM,void *,void *,void *);
389bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *);
3900f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *);
3910f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *);
3927a73cf09SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianAction(DM, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *);
393e3a53471SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, Vec);
3941624f205SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat);
3950f2d7e86SMatthew G. Knepley 
396ef68eab9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *);
397924a1b8fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
398adbe6fbbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *);
399756a1f44SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
400adbe6fbbSMatthew G. Knepley 
4011c41a8caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
4025d16530eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM,Vec,Vec);
403a68b90caSToby Isaac 
404a17985deSToby Isaac /* anchors */
405a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection*, IS*);
406a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS);
407d6a7ad0dSToby Isaac /* tree */
408d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM);
409d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM*);
410dcbd3bf7SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
411da43764aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM*);
412b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt [], PetscInt []);
413b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]);
414d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *);
415d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]);
4166f5f1567SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *);
4173b490a17SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *);
418bbcf1e96SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM,Vec,DM,Vec,PetscSF,PetscSF,PetscInt *,PetscInt *,PetscBool,PetscReal);
419fa534816SMatthew G. Knepley 
420c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMonitorThroughput(DM, void *);
421c0f0dcc3SMatthew G. Knepley 
422fa534816SMatthew G. Knepley /* natural order */
423fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *);
424f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF);
425f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *);
426fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec);
427fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec);
428fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec);
429fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec);
4300bbe5d31Sbarral 
4310bbe5d31Sbarral /* mesh adaptation */
4320f7e110dSbarral PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char [], DM *);
433a8ededdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSnapToGeomModel(DM, PetscInt, const PetscScalar[], PetscScalar[]);
434ca3d3a14SMatthew G. Knepley 
435ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToLocalBasis(DM, Vec);
436ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLocalToGlobalBasis(DM, Vec);
437ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBasisRotation(DM, PetscReal, PetscReal, PetscReal);
438412e9a14SMatthew G. Knepley 
439412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCellRefinerCreate(DM, DMPlexCellRefiner *);
440412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCellRefinerSetUp(DMPlexCellRefiner);
441412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCellRefinerDestroy(DMPlexCellRefiner *);
442412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCellRefinerRefine(DMPlexCellRefiner, DMPolytopeType, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
443412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCellRefinerGetAffineTransforms(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
444412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCellRefinerGetAffineFaceTransforms(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
445412e9a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRefineUniform(DM, DMPlexCellRefiner, DM *);
446552f7358SJed Brown #endif
447