Actual source code: petscds.h
petsc-3.7.3 2016-08-01
1: /*
2: Objects which encapsulate discretizations+continuum residuals
3: */
6: #include <petscfe.h>
7: #include <petscfv.h>
8: #include <petscdstypes.h>
10: PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);
12: PETSC_EXTERN PetscClassId PETSCDS_CLASSID;
14: /*J
15: PetscDSType - String with the name of a PETSc discrete system
17: Level: beginner
19: .seealso: PetscDSSetType(), PetscDS
20: J*/
21: typedef const char *PetscDSType;
22: #define PETSCDSBASIC "basic"
24: typedef void (*PetscPointFunc)(PetscInt, PetscInt, PetscInt,
25: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
26: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
27: PetscReal, const PetscReal[], PetscScalar[]);
28: typedef void (*PetscPointJac)(PetscInt, PetscInt, PetscInt,
29: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
30: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
31: PetscReal, PetscReal, const PetscReal[], PetscScalar[]);
32: typedef void (*PetscBdPointFunc)(PetscInt, PetscInt, PetscInt,
33: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
34: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
35: PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]);
36: typedef void (*PetscBdPointJac)(PetscInt, PetscInt, PetscInt,
37: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
38: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
39: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]);
40: typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *);
43: PETSC_EXTERN PetscFunctionList PetscDSList;
44: PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
45: PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
46: PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
47: PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
48: PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
49: PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
50: PETSC_STATIC_INLINE PetscErrorCode PetscDSViewFromOptions(PetscDS A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
52: PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
53: PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
54: PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
56: PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
57: PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
58: PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
59: PETSC_EXTERN PetscErrorCode PetscDSGetTotalBdDimension(PetscDS, PetscInt *);
60: PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
61: PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
62: PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
63: PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
64: PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
65: PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS, PetscInt *[]);
66: PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
67: PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS, PetscInt *[]);
69: PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
70: PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
71: PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
72: PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
73: PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
74: PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
75: PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
76: PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
77: PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*);
78: PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool, PetscBool);
79: PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
80: void (**)(PetscInt, PetscInt, PetscInt,
81: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
82: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
83: PetscReal, const PetscReal[], PetscScalar[]));
84: PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
85: void (*)(PetscInt, PetscInt, PetscInt,
86: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
87: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
88: PetscReal, const PetscReal[], PetscScalar[]));
89: PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
90: void (**)(PetscInt, PetscInt, PetscInt,
91: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
92: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
93: PetscReal, const PetscReal[], PetscScalar[]),
94: void (**)(PetscInt, PetscInt, PetscInt,
95: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
96: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
97: PetscReal, const PetscReal[], PetscScalar[]));
98: PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
99: void (*)(PetscInt, PetscInt, PetscInt,
100: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
102: PetscReal, const PetscReal[], PetscScalar[]),
103: void (*)(PetscInt, PetscInt, PetscInt,
104: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
105: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
106: PetscReal, const PetscReal[], PetscScalar[]));
107: PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
108: void (**)(PetscInt, PetscInt, PetscInt,
109: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
110: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
111: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
112: void (**)(PetscInt, PetscInt, PetscInt,
113: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
114: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
115: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
116: void (**)(PetscInt, PetscInt, PetscInt,
117: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
118: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
119: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
120: void (**)(PetscInt, PetscInt, PetscInt,
121: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
122: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
123: PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
124: PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
125: void (*)(PetscInt, PetscInt, PetscInt,
126: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
127: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
128: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
129: void (*)(PetscInt, PetscInt, PetscInt,
130: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
131: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
132: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
133: void (*)(PetscInt, PetscInt, PetscInt,
134: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
135: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
136: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
137: void (*)(PetscInt, PetscInt, PetscInt,
138: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
139: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
140: PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
141: PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *);
142: PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
143: void (**)(PetscInt, PetscInt, PetscInt,
144: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
145: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
146: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
147: void (**)(PetscInt, PetscInt, PetscInt,
148: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
149: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
150: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
151: void (**)(PetscInt, PetscInt, PetscInt,
152: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
153: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
154: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
155: void (**)(PetscInt, PetscInt, PetscInt,
156: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
157: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
158: PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
159: PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
160: void (*)(PetscInt, PetscInt, PetscInt,
161: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
162: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
163: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
164: void (*)(PetscInt, PetscInt, PetscInt,
165: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
166: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
167: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
168: void (*)(PetscInt, PetscInt, PetscInt,
169: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
170: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
171: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
172: void (*)(PetscInt, PetscInt, PetscInt,
173: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
174: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
175: PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
176: PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *);
177: PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt,
178: void (**)(PetscInt, PetscInt, PetscInt,
179: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
180: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
181: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
182: void (**)(PetscInt, PetscInt, PetscInt,
183: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
184: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
185: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
186: void (**)(PetscInt, PetscInt, PetscInt,
187: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
188: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
189: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
190: void (**)(PetscInt, PetscInt, PetscInt,
191: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
192: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
193: PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
194: PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt,
195: void (*)(PetscInt, PetscInt, PetscInt,
196: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
197: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
198: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
199: void (*)(PetscInt, PetscInt, PetscInt,
200: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
201: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
202: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
203: void (*)(PetscInt, PetscInt, PetscInt,
204: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
205: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
206: PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
207: void (*)(PetscInt, PetscInt, PetscInt,
208: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
209: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
210: PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
211: PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
212: void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
213: PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
214: void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
215: PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
216: PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
217: PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
218: void (**)(PetscInt, PetscInt, PetscInt,
219: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
220: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
221: PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
222: void (**)(PetscInt, PetscInt, PetscInt,
223: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
224: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
225: PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
226: PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
227: void (*)(PetscInt, PetscInt, PetscInt,
228: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
229: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
230: PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
231: void (*)(PetscInt, PetscInt, PetscInt,
232: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
233: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
234: PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
235: PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
236: void (**)(PetscInt, PetscInt, PetscInt,
237: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
238: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
239: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
240: void (**)(PetscInt, PetscInt, PetscInt,
241: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
242: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
243: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
244: void (**)(PetscInt, PetscInt, PetscInt,
245: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
246: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
247: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
248: void (**)(PetscInt, PetscInt, PetscInt,
249: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
250: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
251: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
252: PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
253: void (*)(PetscInt, PetscInt, PetscInt,
254: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
255: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
256: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
257: void (*)(PetscInt, PetscInt, PetscInt,
258: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
259: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
260: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
261: void (*)(PetscInt, PetscInt, PetscInt,
262: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
263: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
264: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
265: void (*)(PetscInt, PetscInt, PetscInt,
266: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
267: const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
268: PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
269: PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
270: PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
271: PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
272: PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
273: PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
274: PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS);
276: #endif