xref: /petsc/include/petscdmplex.h (revision ee72d991cf2e09b2dca3ed8f78dd34de514b0399)
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 
70c312b8eSJed Brown #include <petscsftypes.h>
8552f7358SJed Brown #include <petscdm.h>
9552f7358SJed Brown 
10552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*);
11830e53efSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, const char[], PetscInt, DM*);
1227c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *);
130fde5641SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*);
1406bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []);
15552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetDimension(DM, PetscInt *);
16552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetDimension(DM, PetscInt);
17552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *);
18552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt);
19552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *);
20552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt);
21552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
22552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
239f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
2477c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
25552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
26552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
27552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
28552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
29552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
30552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
31552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
32552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
33552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
34552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
35552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
36552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
37552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
38d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM);
39552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCoordinateSection(DM, PetscSection *);
40552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCoordinateSection(DM, PetscSection);
41552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetPreallocationCenterDimension(DM,PetscInt);
42cb1e1211SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscSection, PetscSection, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
43552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*);
44552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*);
45552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,const void*);
46552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*);
47552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*);
48552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*);
49552f7358SJed Brown 
50552f7358SJed Brown /*S
51552f7358SJed Brown   DMLabel - Object which encapsulates a subset of the mesh from this DM
52552f7358SJed Brown 
53552f7358SJed Brown   Level: developer
54552f7358SJed Brown 
55552f7358SJed Brown   Concepts: grids, grid refinement
56552f7358SJed Brown 
57552f7358SJed Brown .seealso:  DM, DMPlexCreate(), DMPlexCreateLabel()
58552f7358SJed Brown S*/
59552f7358SJed Brown typedef struct _n_DMLabel *DMLabel;
60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelCreate(const char [], DMLabel *);
61552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelView(DMLabel, PetscViewer);
62552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelDestroy(DMLabel *);
63f807ecf7SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelDuplicate(DMLabel, DMLabel *);
6485de2ee3SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelGetName(DMLabel, const char **);
65552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetValue(DMLabel, PetscInt, PetscInt *);
66552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelSetValue(DMLabel, PetscInt, PetscInt);
67552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelClearValue(DMLabel, PetscInt, PetscInt);
68552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetNumValues(DMLabel, PetscInt *);
69552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetValueIS(DMLabel, IS *);
70552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetStratumSize(DMLabel, PetscInt, PetscInt *);
71552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetStratumIS(DMLabel, PetscInt, IS *);
72552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelClearStratum(DMLabel, PetscInt);
7391d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelCreateIndex(DMLabel, PetscInt, PetscInt);
7491d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelDestroyIndex(DMLabel);
7591d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelHasPoint(DMLabel, PetscInt, PetscBool *);
76552f7358SJed Brown 
77552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateLabel(DM, const char []);
78552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelValue(DM, const char[], PetscInt, PetscInt *);
79552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetLabelValue(DM, const char[], PetscInt, PetscInt);
80552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClearLabelValue(DM, const char[], PetscInt, PetscInt);
81552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelSize(DM, const char[], PetscInt *);
82552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelIdIS(DM, const char[], IS *);
83552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetStratumSize(DM, const char [], PetscInt, PetscInt *);
84552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetStratumIS(DM, const char [], PetscInt, IS *);
85552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClearLabelStratum(DM, const char[], PetscInt);
86552f7358SJed Brown PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection, PetscSF, PetscBool, DMLabel, PetscInt, PetscSection *);
87552f7358SJed Brown 
88552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetNumLabels(DM, PetscInt *);
89552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelName(DM, PetscInt, const char **);
90552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexHasLabel(DM, const char [], PetscBool *);
91552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabel(DM, const char *, DMLabel *);
92552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexAddLabel(DM, DMLabel);
93552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRemoveLabel(DM, const char [], DMLabel *);
94552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
95552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
96552f7358SJed Brown 
97552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
98552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
99552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
100552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
101552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
102552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
103552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
104552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
105552f7358SJed Brown 
10677c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
107ac17c9edSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *);
108552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *);
109552f7358SJed Brown 
110552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *);
111552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
112552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
113552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
114552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
1159db3b87dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []);
1169f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
117930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM);
118552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, const char[], PetscInt, DM*);
119552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexLoad(PetscViewer, DM);
120a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*);
121a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
122e6ccafaeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *);
123552f7358SJed Brown 
124552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []);
125552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, DM *);
126552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm,PetscInt,const PetscInt[],DM *);
127552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
128aa50250dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *);
129552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
130552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
131552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt,const PetscInt [],const PetscInt [], PetscInt,const PetscInt [],const IS [], PetscSection *);
132552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
133bb55d314SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, PetscBool, DM);
13409f723d9SJed Brown PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel);
1352be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
136552f7358SJed Brown 
13718ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */
13818ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
139f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *);
14018ad9376SMatthew G. Knepley 
141834e62ceSMatthew G. Knepley /* FVM Support */
142cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []);
143834e62ceSMatthew G. Knepley 
144552f7358SJed Brown /* FEM Support */
145552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometry(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
1467c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
1477c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
148552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
1497c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
150*ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection);
151552f7358SJed Brown 
152552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
1535bb3eff3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *);
154552f7358SJed Brown 
155552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *);
156a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *);
157552f7358SJed Brown 
158770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
159770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt);
160552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
161552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
162552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
163552f7358SJed Brown 
164552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
165552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
166552f7358SJed Brown 
167552f7358SJed Brown typedef struct {
168552f7358SJed Brown   DM    dm;
169552f7358SJed Brown   Vec   u; /* The base vector for the Jacbobian action J(u) x */
170552f7358SJed Brown   Mat   J; /* Preconditioner for testing */
171552f7358SJed Brown   void *user;
172552f7358SJed Brown } JacActionCtx;
173552f7358SJed Brown 
174a1d24da5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, PetscInt, void (**)(const PetscReal [], PetscScalar *), InsertMode, Vec);
175a1d24da5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFunctionLocal(DM, PetscInt, void (**)(const PetscReal [], PetscScalar *), InsertMode, Vec);
176a1d24da5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2Diff(DM, PetscQuadrature[], void (**)(const PetscReal [], PetscScalar *), Vec, PetscReal *);
177552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeResidualFEM(DM, Vec, Vec, void *);
178552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianActionFEM(DM, Mat, Vec, Vec, void *);
179552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianFEM(DM, Vec, Mat, Mat, MatStructure*,void *);
180552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetFEMIntegration(DM,
181552f7358SJed Brown                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
182552f7358SJed Brown                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
183552f7358SJed Brown                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
184552f7358SJed Brown                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]),
18502a80efeSMatthew G. Knepley                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
18602a80efeSMatthew G. Knepley                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
18702a80efeSMatthew G. Knepley                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
18802a80efeSMatthew G. Knepley                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), PetscScalar[]),
189552f7358SJed Brown                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[],
190552f7358SJed Brown                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
191552f7358SJed Brown                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
192552f7358SJed Brown                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
193552f7358SJed Brown                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
194552f7358SJed Brown                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]),
195552f7358SJed Brown                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
196552f7358SJed Brown                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
197552f7358SJed Brown                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
198552f7358SJed Brown                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
199552f7358SJed Brown                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
200552f7358SJed Brown                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]));
201552f7358SJed Brown #endif
202