xref: /petsc/include/petscdmplex.h (revision 03da9461e5c4ceb1d4e5e8877d1e025447eb954a)
1552f7358SJed Brown /*
2552f7358SJed Brown   DMPlex, for parallel unstructured distributed mesh problems.
3552f7358SJed Brown */
4552f7358SJed Brown #if !defined(__PETSCDMPLEX_H)
5552f7358SJed Brown #define __PETSCDMPLEX_H
6552f7358SJed Brown 
7552f7358SJed Brown #include <petscdm.h>
8a0845e3aSMatthew G. Knepley #include <petscdt.h>
9a0845e3aSMatthew G. Knepley #include <petscfe.h>
100bacfa0aSMatthew G. Knepley #include <petscfv.h>
11799f573fSMatthew G. Knepley #include <petscsftypes.h>
12c330f8ffSToby Isaac #include <petscdmfield.h>
13552f7358SJed Brown 
1477623264SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCPARTITIONER_CLASSID;
1577623264SMatthew G. Knepley 
1677623264SMatthew G. Knepley /*J
1777623264SMatthew G. Knepley   PetscPartitionerType - String with the name of a PETSc graph partitioner
1877623264SMatthew G. Knepley 
1977623264SMatthew G. Knepley   Level: beginner
2077623264SMatthew G. Knepley 
2177623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitioner
2277623264SMatthew G. Knepley J*/
2377623264SMatthew G. Knepley typedef const char *PetscPartitionerType;
2477623264SMatthew G. Knepley #define PETSCPARTITIONERCHACO    "chaco"
2577623264SMatthew G. Knepley #define PETSCPARTITIONERPARMETIS "parmetis"
26137cd93aSLisandro Dalcin #define PETSCPARTITIONERPTSCOTCH "ptscotch"
2777623264SMatthew G. Knepley #define PETSCPARTITIONERSHELL    "shell"
28555a9cf8SMatthew G. Knepley #define PETSCPARTITIONERSIMPLE   "simple"
296462276dSToby Isaac #define PETSCPARTITIONERGATHER   "gather"
30de68236aSVaclav Hapla #define PETSCPARTITIONERMATPARTITIONING "matpartitioning"
3177623264SMatthew G. Knepley 
3277623264SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscPartitionerList;
3377623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate(MPI_Comm, PetscPartitioner *);
3477623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *);
3577623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetType(PetscPartitioner, PetscPartitionerType);
3677623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerGetType(PetscPartitioner, PetscPartitionerType *);
3777623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetUp(PetscPartitioner);
3877623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner);
398aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
4077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerView(PetscPartitioner, PetscViewer);
4177623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerRegister(const char [], PetscErrorCode (*)(PetscPartitioner));
4277623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterDestroy(void);
4377623264SMatthew G. Knepley 
44f8987ae8SMichael Lange PETSC_EXTERN PetscErrorCode PetscPartitionerPartition(PetscPartitioner, DM, PetscSection, IS *);
4577623264SMatthew G. Knepley 
465680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner, PetscInt, const PetscInt[], const PetscInt[]);
47df8115dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner, PetscBool);
48df8115dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner, PetscBool *);
4977623264SMatthew G. Knepley 
50de68236aSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning(PetscPartitioner part, MatPartitioning *mp);
51de68236aSVaclav Hapla 
52552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*);
5327c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *);
540fde5641SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*);
55fba955ccSBarry Smith PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const PetscReal[], PetscSF *, DM *);
5606bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []);
578415267dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, PetscInt, PetscBool, DM*);
58a9074c1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetOptionsPrefix(DM, const char []);
59552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *);
60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt);
61552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *);
62552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt);
63feb68dd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt);
64552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
65552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
669f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
6777c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
68552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
69552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
70552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
71552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
72552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
73552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
74552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
75552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
768cb4d582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *);
77552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
78552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
79552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
80552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
81552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
824e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*);
83cb11e54dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexReverseCell(DM, PetscInt);
84d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM);
8575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
864e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *);
87e101f074SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
88552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*);
89081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,void*);
90a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*);
91a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*);
92a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*);
93081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,void*);
94552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*);
95552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*);
96a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*);
97a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*);
98a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*);
99081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*, void*);
100552f7358SJed Brown 
1018b51c812SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *);
102552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
103552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
104ef48cebcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *);
10508a22f4bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *);
106552f7358SJed Brown 
10775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
1083ded2ed9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *);
10975d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
11075d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
11175d3a19aSMatthew G. Knepley 
11275d3a19aSMatthew G. Knepley /* Topological Operations */
113552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
114552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
115552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
116552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
117552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
118552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
119552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
120552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
121552f7358SJed Brown 
12275d3a19aSMatthew G. Knepley /* Mesh Generation */
123552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *);
1243a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegister(const char[],PetscErrorCode (*)(DM,PetscBool,DM*),PetscErrorCode (*)(DM,double*,DM*),PetscInt);
1253a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegisterAll(void);
1263a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegisterDestroy(void);
127930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM);
1281df5d5c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscBool, PetscReal, DM *);
12926492d91SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSquareBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []);
130552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []);
131768d5fceSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
13265a81367SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, DM *);
133dbc1dc17SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, PetscInt, DMBoundaryType, DM *);
13424119c2aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *);
13500dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *);
13600dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, DM *);
137552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
138459e96c1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []);
139ca8062c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM);
14058723a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscBool, PetscInt);
1419bf0dad6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscBool, PetscInt);
142*03da9461SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM);
14375d3a19aSMatthew G. Knepley 
144d9deefdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *);
145776fd405SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *);
146d9deefdfSMatthew G. Knepley 
147d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], PetscBool, DM *);
148d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
149d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char [], PetscBool, DM *);
150d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *);
151d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *);
152d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *);
153d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char [], PetscBool, DM *);
154d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *);
155d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char [], PetscBool, DM *);
156d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char [], PetscBool, DM *);
157f2801cd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char [], PetscBool, DM *);
158d6218766SMatthew G. Knepley 
15975d3a19aSMatthew G. Knepley /* Mesh Partitioning and Distribution */
160552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
1615680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *);
16271bb2955SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner);
163ac17c9edSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *);
1643fa7752dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *);
165552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *);
1661fd9873aSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel);
1675abbe4feSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel);
16824d039d7SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel);
169be200f8dSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel);
170aa3148a8SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *);
17198ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool);
17298ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *);
17380cf41d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF*, DM*);
174b9f40539SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *);
175552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM,PetscSF,PetscSection,Vec,PetscSection,Vec);
176b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *);
177552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM,PetscSF,PetscSection,MPI_Datatype,void*,PetscSection,void**);
17845800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM);
1796462276dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, DM*);
1803c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM,PetscInt,PetscInt*,PetscInt[],void*),void*);
1813c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM,PetscInt,PetscInt*,PetscInt[],void*),void**);
18270034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseCone(DM,PetscBool);
18370034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseCone(DM,PetscBool*);
18470034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseClosure(DM,PetscBool);
18570034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseClosure(DM,PetscBool*);
186a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM,PetscBool);
187a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM,PetscBool*);
18870034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]);
189f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF);
190f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *);
19175d3a19aSMatthew G. Knepley 
192c99f2573SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *);
193f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *);
194f5cedc29SMatthew G. Knepley 
195b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *);
196b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *);
197b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *);
198a9f1d5b2SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlap(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *);
19946f9b1c3SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *);
20045800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *);
201b0a623aaSMatthew G. Knepley 
20275d3a19aSMatthew G. Knepley /* Submesh Support */
203158acfadSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM*);
2047db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, DMLabel *, DMLabel *, DM *, DM *);
205552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*);
206552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
207552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *);
208559a1558SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpoint(DM, PetscInt, PetscInt *);
209552f7358SJed Brown 
2102be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
2110f66a230SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscBool, DM);
2126cf0e42fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel);
213f402d5e4SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel);
214552f7358SJed Brown 
21575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
21675d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
21775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
21875d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
2191e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal [], PetscReal *));
2201e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal [], PetscReal *));
2212389894bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *);
2220aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *);
2230aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool);
224e5337592SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexRefineSimplexToTensor(DM, DM*);
225552f7358SJed Brown 
22618ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */
22718ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
228f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *);
22918ad9376SMatthew G. Knepley 
2300f0bf202SMatthew G. Knepley /* Geometry Support */
2310f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *);
2320f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal);
233741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]);
234741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]);
235741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]);
2360f0bf202SMatthew G. Knepley 
237c4eade1cSMatthew G. Knepley /* Point Location */
238c4eade1cSMatthew G. Knepley typedef struct _PetscGridHash *PetscGridHash;
239c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *);
240c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar []);
241c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt [], const PetscReal []);
2421c6dfc3eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar [], PetscInt [], PetscInt []);
243c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *);
244c4eade1cSMatthew G. Knepley 
245834e62ceSMatthew G. Knepley /* FVM Support */
246cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []);
2470f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *);
248856ac710SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *);
249b27d5b9eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *);
250834e62ceSMatthew G. Knepley 
251552f7358SJed Brown /* FEM Support */
252cf0b7c11SKarl Rupp PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *);
253bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
2541c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[],
255b278463cSMatthew G. Knepley                                                                 PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec );
2561c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
257b278463cSMatthew G. Knepley                                                                      void (*)(PetscInt, PetscInt, PetscInt,
258b278463cSMatthew G. Knepley                                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
259b278463cSMatthew G. Knepley                                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
260191494d9SMatthew G. Knepley                                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[],
261b278463cSMatthew G. Knepley                                                                               PetscScalar[]),
262b278463cSMatthew G. Knepley                                                                      void *, Vec);
2631c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM, PetscReal, Vec, Vec, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [],
264b278463cSMatthew G. Knepley                                                               PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *, Vec);
265e752be1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel);
266d7ddef95SMatthew G. Knepley 
267ab1d0545SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt, const PetscInt [], const PetscInt [], PetscInt, const PetscInt [], const IS [], const IS [], IS, PetscSection *);
268be36d101SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection*);
26975d3a19aSMatthew G. Knepley 
2708e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
271dfccc68fSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
272d6143a4eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]);
2739d150b73SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]);
274d6143a4eSToby Isaac 
2757c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
2767c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
277552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
2787c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
2796ecaa68aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscInt *);
28046bdb399SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscInt *, PetscInt **,PetscInt *);
2810405ed22SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
2827c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]);
283ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection);
2847391a63aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSpectralClosurePermutation(DM, PetscInt, PetscSection);
285552f7358SJed Brown 
286552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *);
2877db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *);
288552f7358SJed Brown 
289770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
290770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt);
291552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
292552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
293552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
294552f7358SJed Brown 
295552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
296552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
297552f7358SJed Brown 
298552f7358SJed Brown typedef struct {
299552f7358SJed Brown   DM    dm;
300552f7358SJed Brown   Vec   u; /* The base vector for the Jacbobian action J(u) x */
301552f7358SJed Brown   Mat   J; /* Preconditioner for testing */
302552f7358SJed Brown   void *user;
303552f7358SJed Brown } JacActionCtx;
304552f7358SJed Brown 
3053351dd3dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM, Vec);
306b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt);
307b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt*);
308e29c0834SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec,
309e29c0834SMatthew G. Knepley                                                     void (**)(PetscInt, PetscInt, PetscInt,
310e29c0834SMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
311e29c0834SMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
312e29c0834SMatthew G. Knepley                                                               PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec);
313574a98acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
3140163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]);
3150163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, Vec);
316338f77d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *);
317b8feb594SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *);
3189b6f715bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdIntegral(DM, Vec, DMLabel, PetscInt, const PetscInt[],
3199b6f715bSMatthew G. Knepley                                                     void (*)(PetscInt, PetscInt, PetscInt,
3209b6f715bSMatthew G. Knepley                                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3219b6f715bSMatthew G. Knepley                                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3229b6f715bSMatthew G. Knepley                                                              PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
3239b6f715bSMatthew G. Knepley                                                     PetscScalar *, void *);
32468132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, Mat, void *);
32568132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *);
3261555c271SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec);
3277c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *);
3282cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *);
3292cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *);
3301c41a8caSMatthew G. Knepley 
331026175e5SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, MatNullSpace *);
332fd86e548SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *);
333026175e5SToby Isaac 
3349f520fc2SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM,void *,void *,void *);
335bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *);
3360f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *);
3370f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *);
3387a73cf09SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianAction(DM, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *);
339e3a53471SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, Vec);
3401624f205SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat);
3410f2d7e86SMatthew G. Knepley 
342ef68eab9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *);
343924a1b8fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
344adbe6fbbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *);
345756a1f44SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
346adbe6fbbSMatthew G. Knepley 
3471c41a8caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
3485d16530eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM,Vec,Vec);
349a68b90caSToby Isaac 
350a17985deSToby Isaac /* anchors */
351a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection*, IS*);
352a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS);
353d6a7ad0dSToby Isaac /* tree */
354d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM);
355d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM*);
356dcbd3bf7SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
357da43764aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM*);
358b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt [], PetscInt []);
359b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]);
360d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *);
361d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]);
3626f5f1567SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *);
3633b490a17SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *);
364bbcf1e96SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM,Vec,DM,Vec,PetscSF,PetscSF,PetscInt *,PetscInt *,PetscBool,PetscReal);
365fa534816SMatthew G. Knepley 
366fa534816SMatthew G. Knepley /* natural order */
367fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *);
368f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF);
369f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *);
370fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec);
371fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec);
372fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec);
373fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec);
3740bbe5d31Sbarral 
3750bbe5d31Sbarral /* mesh adaptation */
3760f7e110dSbarral PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char [], DM *);
377552f7358SJed Brown #endif
378