Actual source code: petscds.h

  1: /*
  2:       Objects which encapsulate discretizations+continuum residuals
  3: */
  4: #if !defined(PETSCDS_H)
  5: #define PETSCDS_H
  6: #include <petscfe.h>
  7: #include <petscfv.h>
  8: #include <petscdstypes.h>

 10: PETSC_EXTERN PetscClassId PETSCWEAKFORM_CLASSID;

 12: PETSC_EXTERN PetscErrorCode PetscWeakFormCreate(MPI_Comm, PetscWeakForm *);
 13: PETSC_EXTERN PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *);
 14: PETSC_EXTERN PetscErrorCode PetscWeakFormView(PetscWeakForm, PetscViewer);
 15: PETSC_EXTERN PetscErrorCode PetscWeakFormCopy(PetscWeakForm, PetscWeakForm);
 16: PETSC_EXTERN PetscErrorCode PetscWeakFormClear(PetscWeakForm);
 17: PETSC_EXTERN PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm, PetscInt *);
 18: PETSC_EXTERN PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm, PetscInt);
 19: PETSC_EXTERN PetscErrorCode PetscFormKeySort(PetscInt, PetscFormKey[]);
 20: PETSC_EXTERN PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm, DMLabel, PetscInt, const PetscInt[]);
 21: PETSC_EXTERN PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm, DMLabel);
 22: PETSC_EXTERN PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscWeakFormKind, PetscInt);

 24: PETSC_EXTERN PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt *,
 25:                                                       void (***)(PetscInt, PetscInt, PetscInt,
 26:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 27:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 28:                                                                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 29: PETSC_EXTERN PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
 30:                                                       void (*)(PetscInt, PetscInt, PetscInt,
 31:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 32:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 33:                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 34: PETSC_EXTERN PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
 35:                                                       void (**)(PetscInt, PetscInt, PetscInt,
 36:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 37:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 38:                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 39: PETSC_EXTERN PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
 40:                                                       void (**)(PetscInt, PetscInt, PetscInt,
 41:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 42:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 43:                                                                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 44: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
 45:                                                       void (*)(PetscInt, PetscInt, PetscInt,
 46:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 47:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 48:                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 49: PETSC_EXTERN PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
 50:                                                      PetscInt *,
 51:                                                      void (***)(PetscInt, PetscInt, PetscInt,
 52:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 53:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 54:                                                                 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
 55:                                                      PetscInt *,
 56:                                                      void (***)(PetscInt, PetscInt, PetscInt,
 57:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 58:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 59:                                                                 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 60: PETSC_EXTERN PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
 61:                                                      void (*)(PetscInt, PetscInt, PetscInt,
 62:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 63:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 64:                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
 65:                                                       void (*)(PetscInt, PetscInt, PetscInt,
 66:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 67:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 68:                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 69: PETSC_EXTERN PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
 70:                                                      PetscInt,
 71:                                                      void (**)(PetscInt, PetscInt, PetscInt,
 72:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 73:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 74:                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
 75:                                                      PetscInt,
 76:                                                      void (**)(PetscInt, PetscInt, PetscInt,
 77:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 78:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 79:                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 80: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
 81:                                                      PetscInt,
 82:                                                      void (*)(PetscInt, PetscInt, PetscInt,
 83:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 84:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 85:                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
 86:                                                      PetscInt,
 87:                                                      void (*)(PetscInt, PetscInt, PetscInt,
 88:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 89:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 90:                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
 91: PETSC_EXTERN PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm, PetscBool *);
 92: PETSC_EXTERN PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
 93:                                                      PetscInt *,
 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, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
 98:                                                      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, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
103:                                                      PetscInt *,
104:                                                      void (***)(PetscInt, PetscInt, PetscInt,
105:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
106:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
107:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
108:                                                      PetscInt *,
109:                                                      void (***)(PetscInt, PetscInt, PetscInt,
110:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
111:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
112:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
113: PETSC_EXTERN PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
114:                                                      void (*)(PetscInt, PetscInt, PetscInt,
115:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
116:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
117:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
118:                                                       void (*)(PetscInt, PetscInt, PetscInt,
119:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
120:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
121:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
122:                                                       void (*)(PetscInt, PetscInt, PetscInt,
123:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
124:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
125:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
126:                                                       void (*)(PetscInt, PetscInt, PetscInt,
127:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
128:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
129:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
130: PETSC_EXTERN PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
131:                                                      PetscInt,
132:                                                      void (**)(PetscInt, PetscInt, PetscInt,
133:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
134:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
135:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
136:                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
141:                                                      PetscInt,
142:                                                      void (**)(PetscInt, PetscInt, PetscInt,
143:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
144:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
145:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
146:                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]));
151: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
152:                                                      PetscInt,
153:                                                      void (*)(PetscInt, PetscInt, PetscInt,
154:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
155:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
156:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
157:                                                      PetscInt,
158:                                                      void (*)(PetscInt, PetscInt, PetscInt,
159:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
160:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
161:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
162:                                                      PetscInt,
163:                                                      void (*)(PetscInt, PetscInt, PetscInt,
164:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
165:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
166:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
167:                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]));
172: PETSC_EXTERN PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm, PetscBool *);
173: PETSC_EXTERN PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
174:                                                      PetscInt *,
175:                                                      void (***)(PetscInt, PetscInt, PetscInt,
176:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
177:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
178:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
179:                                                      PetscInt *,
180:                                                      void (***)(PetscInt, PetscInt, PetscInt,
181:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
182:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
183:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
184:                                                      PetscInt *,
185:                                                      void (***)(PetscInt, PetscInt, PetscInt,
186:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
187:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
188:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
189:                                                      PetscInt *,
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[], PetscInt, const PetscScalar[], PetscScalar[]));
194: PETSC_EXTERN PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, 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[], PetscInt, const PetscScalar[], 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[], PetscInt, const PetscScalar[], 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[], PetscInt, const PetscScalar[], 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[], PetscInt, const PetscScalar[], PetscScalar[]));
211: PETSC_EXTERN PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
212:                                                      PetscInt,
213:                                                      void (**)(PetscInt, PetscInt, PetscInt,
214:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
215:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
216:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
217:                                                      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, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
222:                                                      PetscInt,
223:                                                      void (**)(PetscInt, PetscInt, PetscInt,
224:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
225:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
226:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
227:                                                      PetscInt,
228:                                                      void (**)(PetscInt, PetscInt, PetscInt,
229:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
230:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
231:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
232: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
233:                                                      PetscInt,
234:                                                      void (*)(PetscInt, PetscInt, PetscInt,
235:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
236:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
237:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
238:                                                      PetscInt,
239:                                                      void (*)(PetscInt, PetscInt, PetscInt,
240:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
241:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
242:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
243:                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
248:                                                      PetscInt,
249:                                                      void (*)(PetscInt, PetscInt, PetscInt,
250:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
251:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
252:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
253: PETSC_EXTERN PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm, PetscBool *);
254: PETSC_EXTERN PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
255:                                                      PetscInt *,
256:                                                      void (***)(PetscInt, PetscInt, PetscInt,
257:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
258:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
259:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
260:                                                      PetscInt *,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
265:                                                      PetscInt *,
266:                                                      void (***)(PetscInt, PetscInt, PetscInt,
267:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
268:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
269:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
270:                                                      PetscInt *,
271:                                                      void (***)(PetscInt, PetscInt, PetscInt,
272:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
273:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
274:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
275: PETSC_EXTERN PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
276:                                                      void (*)(PetscInt, PetscInt, PetscInt,
277:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
278:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
279:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
280:                                                       void (*)(PetscInt, PetscInt, PetscInt,
281:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
282:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
283:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
284:                                                       void (*)(PetscInt, PetscInt, PetscInt,
285:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
286:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
287:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
288:                                                       void (*)(PetscInt, PetscInt, PetscInt,
289:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
290:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
291:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
292: PETSC_EXTERN PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
293:                                                      PetscInt,
294:                                                      void (**)(PetscInt, PetscInt, PetscInt,
295:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
296:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
297:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
298:                                                      PetscInt,
299:                                                      void (**)(PetscInt, PetscInt, PetscInt,
300:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
301:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
302:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
303:                                                      PetscInt,
304:                                                      void (**)(PetscInt, PetscInt, PetscInt,
305:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
306:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
307:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
308:                                                      PetscInt,
309:                                                      void (**)(PetscInt, PetscInt, PetscInt,
310:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
311:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
312:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
313: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
314:                                                      PetscInt,
315:                                                      void (*)(PetscInt, PetscInt, PetscInt,
316:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
317:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
318:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
319:                                                      PetscInt,
320:                                                      void (*)(PetscInt, PetscInt, PetscInt,
321:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
322:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
323:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
324:                                                      PetscInt,
325:                                                      void (*)(PetscInt, PetscInt, PetscInt,
326:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
327:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
328:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
329:                                                      PetscInt,
330:                                                      void (*)(PetscInt, PetscInt, PetscInt,
331:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
332:                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
333:                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
334: PETSC_EXTERN PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
335:                                                        PetscInt *,
336:                                                        void (***)(PetscInt, PetscInt, PetscInt,
337:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
338:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
339:                                                                   PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
340:                                                        PetscInt *,
341:                                                        void (***)(PetscInt, PetscInt, PetscInt,
342:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
343:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
344:                                                                   PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
345: PETSC_EXTERN PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
346:                                                        void (*)(PetscInt, PetscInt, PetscInt,
347:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
348:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
349:                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
350:                                                        void (*)(PetscInt, PetscInt, PetscInt,
351:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
352:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
353:                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
354: PETSC_EXTERN PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
355:                                                        PetscInt,
356:                                                        void (**)(PetscInt, PetscInt, PetscInt,
357:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
358:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
359:                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
360:                                                        PetscInt,
361:                                                        void (**)(PetscInt, PetscInt, PetscInt,
362:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
363:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
364:                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
365: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
366:                                                        PetscInt,
367:                                                        void (*)(PetscInt, PetscInt, PetscInt,
368:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
369:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
370:                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
371:                                                        PetscInt,
372:                                                        void (*)(PetscInt, PetscInt, PetscInt,
373:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
374:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
375:                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
376: PETSC_EXTERN PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm, PetscBool *);
377: PETSC_EXTERN PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
378:                                                        PetscInt *,
379:                                                        void (***)(PetscInt, PetscInt, PetscInt,
380:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
381:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
382:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
383:                                                        PetscInt *,
384:                                                        void (***)(PetscInt, PetscInt, PetscInt,
385:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
386:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
387:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
388:                                                        PetscInt *,
389:                                                        void (***)(PetscInt, PetscInt, PetscInt,
390:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
391:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
392:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
393:                                                        PetscInt *,
394:                                                        void (***)(PetscInt, PetscInt, PetscInt,
395:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
396:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
397:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
398: PETSC_EXTERN PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
399:                                                        void (*)(PetscInt, PetscInt, PetscInt,
400:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
401:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
402:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
403:                                                        void (*)(PetscInt, PetscInt, PetscInt,
404:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
405:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
406:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
407:                                                        void (*)(PetscInt, PetscInt, PetscInt,
408:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
409:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
410:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
411:                                                        void (*)(PetscInt, PetscInt, PetscInt,
412:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
413:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
414:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
415: PETSC_EXTERN PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
416:                                                        PetscInt,
417:                                                        void (**)(PetscInt, PetscInt, PetscInt,
418:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
419:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
420:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
421:                                                        PetscInt,
422:                                                        void (**)(PetscInt, PetscInt, PetscInt,
423:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
424:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
425:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
426:                                                        PetscInt,
427:                                                        void (**)(PetscInt, PetscInt, PetscInt,
428:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
429:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
430:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
431:                                                        PetscInt,
432:                                                        void (**)(PetscInt, PetscInt, PetscInt,
433:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
434:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
435:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
436: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
437:                                                        PetscInt,
438:                                                        void (*)(PetscInt, PetscInt, PetscInt,
439:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
440:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
441:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
442:                                                        PetscInt,
443:                                                        void (*)(PetscInt, PetscInt, PetscInt,
444:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
445:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
446:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
447:                                                        PetscInt,
448:                                                        void (*)(PetscInt, PetscInt, PetscInt,
449:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
450:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
451:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
452:                                                        PetscInt,
453:                                                        void (*)(PetscInt, PetscInt, PetscInt,
454:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
455:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
456:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
457: PETSC_EXTERN PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm, PetscBool *);
458: PETSC_EXTERN PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
459:                                                        PetscInt *,
460:                                                        void (***)(PetscInt, PetscInt, PetscInt,
461:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
462:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
463:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
464:                                                        PetscInt *,
465:                                                        void (***)(PetscInt, PetscInt, PetscInt,
466:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
467:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
468:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
469:                                                        PetscInt *,
470:                                                        void (***)(PetscInt, PetscInt, PetscInt,
471:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
472:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
473:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
474:                                                        PetscInt *,
475:                                                        void (***)(PetscInt, PetscInt, PetscInt,
476:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
477:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
478:                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
479: PETSC_EXTERN PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
480:                                                        void (*)(PetscInt, PetscInt, PetscInt,
481:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
482:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
483:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
484:                                                        void (*)(PetscInt, PetscInt, PetscInt,
485:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
486:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
487:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
488:                                                        void (*)(PetscInt, PetscInt, PetscInt,
489:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
490:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
491:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
492:                                                        void (*)(PetscInt, PetscInt, PetscInt,
493:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
494:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
495:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
496: PETSC_EXTERN PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
497:                                                        PetscInt,
498:                                                        void (**)(PetscInt, PetscInt, PetscInt,
499:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
500:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
501:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
502:                                                        PetscInt,
503:                                                        void (**)(PetscInt, PetscInt, PetscInt,
504:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
505:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
506:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
507:                                                        PetscInt,
508:                                                        void (**)(PetscInt, PetscInt, PetscInt,
509:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
510:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
511:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
512:                                                        PetscInt,
513:                                                        void (**)(PetscInt, PetscInt, PetscInt,
514:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
515:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
516:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
517: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
518:                                                        PetscInt,
519:                                                        void (*)(PetscInt, PetscInt, PetscInt,
520:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
521:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
522:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
523:                                                        PetscInt,
524:                                                        void (*)(PetscInt, PetscInt, PetscInt,
525:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
526:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
527:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
528:                                                        PetscInt,
529:                                                        void (*)(PetscInt, PetscInt, PetscInt,
530:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
531:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
532:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
533:                                                        PetscInt,
534:                                                        void (*)(PetscInt, PetscInt, PetscInt,
535:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
536:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
537:                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
538: PETSC_EXTERN PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
539:                                                           PetscInt *,
540:                                                           void (***)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
541: PETSC_EXTERN PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
542:                                                           PetscInt,
543:                                                           void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
544: PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
545:                                                           PetscInt,
546:                                                           void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));

548: PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);

550: PETSC_EXTERN PetscClassId PETSCDS_CLASSID;

552: /*J
553:   PetscDSType - String with the name of a PETSc discrete system

555:   Level: beginner

557: .seealso: PetscDSSetType(), PetscDS
558: J*/
559: typedef const char *PetscDSType;
560: #define PETSCDSBASIC "basic"

562: typedef enum {PETSC_DISC_NONE, PETSC_DISC_FE, PETSC_DISC_FV} PetscDiscType;

564: typedef void (*PetscPointFunc)(PetscInt, PetscInt, PetscInt,
565:                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
566:                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
567:                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
568: typedef void (*PetscPointJac)(PetscInt, PetscInt, PetscInt,
569:                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
570:                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
571:                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
572: typedef void (*PetscBdPointFunc)(PetscInt, PetscInt, PetscInt,
573:                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
574:                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
575:                                  PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
576: typedef void (*PetscBdPointJac)(PetscInt, PetscInt, PetscInt,
577:                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
578:                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
579:                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
580: typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
581: typedef PetscErrorCode (*PetscSimplePointFunc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);

583: PETSC_EXTERN PetscFunctionList PetscDSList;
584: PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
585: PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
586: PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
587: PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
588: PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
589: PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
590: PETSC_EXTERN PetscErrorCode PetscDSViewFromOptions(PetscDS,PetscObject,const char[]);

592: PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
593: PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
594: PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);

596: PETSC_EXTERN PetscErrorCode PetscDSGetHeightSubspace(PetscDS, PetscInt, PetscDS *);
597: PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
598: PETSC_EXTERN PetscErrorCode PetscDSGetCoordinateDimension(PetscDS, PetscInt *);
599: PETSC_EXTERN PetscErrorCode PetscDSSetCoordinateDimension(PetscDS, PetscInt);
600: PETSC_EXTERN PetscErrorCode PetscDSIsCohesive(PetscDS, PetscBool *);
601: PETSC_EXTERN PetscErrorCode PetscDSGetNumCohesive(PetscDS, PetscInt *);
602: PETSC_EXTERN PetscErrorCode PetscDSGetCohesive(PetscDS, PetscInt, PetscBool *);
603: PETSC_EXTERN PetscErrorCode PetscDSSetCohesive(PetscDS, PetscInt, PetscBool);
604: PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
605: PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
606: PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
607: PETSC_EXTERN PetscErrorCode PetscDSGetFieldIndex(PetscDS, PetscObject, PetscInt *);
608: PETSC_EXTERN PetscErrorCode PetscDSGetFieldSize(PetscDS, PetscInt, PetscInt *);
609: PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
610: PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffsetCohesive(PetscDS, PetscInt, PetscInt *);
611: PETSC_EXTERN PetscErrorCode PetscDSGetDimensions(PetscDS, PetscInt *[]);
612: PETSC_EXTERN PetscErrorCode PetscDSGetComponents(PetscDS, PetscInt *[]);
613: PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
614: PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
615: PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
616: PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsetsCohesive(PetscDS, PetscInt, PetscInt *[]);
617: PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsetsCohesive(PetscDS, PetscInt, PetscInt *[]);

619: PETSC_EXTERN PetscErrorCode PetscDSGetWeakForm(PetscDS, PetscWeakForm *);
620: PETSC_EXTERN PetscErrorCode PetscDSSetWeakForm(PetscDS, PetscWeakForm);
621: PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
622: PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
623: PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
624: PETSC_EXTERN PetscErrorCode PetscDSGetQuadrature(PetscDS, PetscQuadrature*);
625: PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
626: PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
627: PETSC_EXTERN PetscErrorCode PetscDSGetJetDegree(PetscDS, PetscInt, PetscInt*);
628: PETSC_EXTERN PetscErrorCode PetscDSSetJetDegree(PetscDS, PetscInt, PetscInt);
629: PETSC_EXTERN PetscErrorCode PetscDSGetConstants(PetscDS, PetscInt *, const PetscScalar *[]);
630: PETSC_EXTERN PetscErrorCode PetscDSSetConstants(PetscDS, PetscInt, PetscScalar[]);
631: PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
632:                                                 void (**)(PetscInt, PetscInt, PetscInt,
633:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
634:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
635:                                                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
636: PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
637:                                                 void (*)(PetscInt, PetscInt, PetscInt,
638:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
639:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
640:                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
641: PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
642:                                                void (**)(PetscInt, PetscInt, PetscInt,
643:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
644:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
645:                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
646:                                                void (**)(PetscInt, PetscInt, PetscInt,
647:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
648:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
649:                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
650: PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
651:                                                void (*)(PetscInt, PetscInt, PetscInt,
652:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
653:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
654:                                                         PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
655:                                                void (*)(PetscInt, PetscInt, PetscInt,
656:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
657:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
658:                                                         PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
659: PETSC_EXTERN PetscErrorCode PetscDSGetRHSResidual(PetscDS, PetscInt,
660:                                                   void (**)(PetscInt, PetscInt, PetscInt,
661:                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
662:                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
663:                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
664:                                                   void (**)(PetscInt, PetscInt, PetscInt,
665:                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
666:                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
667:                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
668: PETSC_EXTERN PetscErrorCode PetscDSSetRHSResidual(PetscDS, PetscInt,
669:                                                   void (*)(PetscInt, PetscInt, PetscInt,
670:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
671:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
672:                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
673:                                                   void (*)(PetscInt, PetscInt, PetscInt,
674:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
675:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
676:                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
677: PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *);
678: PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
679:                                                void (**)(PetscInt, PetscInt, PetscInt,
680:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
681:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
682:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
683:                                                void (**)(PetscInt, PetscInt, PetscInt,
684:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
685:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
686:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
687:                                                void (**)(PetscInt, PetscInt, PetscInt,
688:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
689:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
690:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
691:                                                void (**)(PetscInt, PetscInt, PetscInt,
692:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
693:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
694:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
695: PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
696:                                                void (*)(PetscInt, PetscInt, PetscInt,
697:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
698:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
699:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
700:                                                void (*)(PetscInt, PetscInt, PetscInt,
701:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
702:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
703:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
704:                                                void (*)(PetscInt, PetscInt, PetscInt,
705:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
706:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
707:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
708:                                                void (*)(PetscInt, PetscInt, PetscInt,
709:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
710:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
711:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
712: PETSC_EXTERN PetscErrorCode PetscDSUseJacobianPreconditioner(PetscDS, PetscBool);
713: PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *);
714: PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
715:                                                void (**)(PetscInt, PetscInt, PetscInt,
716:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
717:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
718:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
719:                                                void (**)(PetscInt, PetscInt, PetscInt,
720:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
721:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
722:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
723:                                                void (**)(PetscInt, PetscInt, PetscInt,
724:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
725:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
726:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
727:                                                void (**)(PetscInt, PetscInt, PetscInt,
728:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
729:                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
730:                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
731: PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
732:                                                void (*)(PetscInt, PetscInt, PetscInt,
733:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
734:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
735:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
736:                                                void (*)(PetscInt, PetscInt, PetscInt,
737:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
738:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
739:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
740:                                                void (*)(PetscInt, PetscInt, PetscInt,
741:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
742:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
743:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
744:                                                void (*)(PetscInt, PetscInt, PetscInt,
745:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
746:                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
747:                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
748: PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *);
749: PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt,
750:                                                       void (**)(PetscInt, PetscInt, PetscInt,
751:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
752:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
753:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
754:                                                       void (**)(PetscInt, PetscInt, PetscInt,
755:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
756:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
757:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
758:                                                       void (**)(PetscInt, PetscInt, PetscInt,
759:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
760:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
761:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
762:                                                       void (**)(PetscInt, PetscInt, PetscInt,
763:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
764:                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
765:                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
766: PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt,
767:                                                       void (*)(PetscInt, PetscInt, PetscInt,
768:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
769:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
770:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
771:                                                       void (*)(PetscInt, PetscInt, PetscInt,
772:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
773:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
774:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
775:                                                       void (*)(PetscInt, PetscInt, PetscInt,
776:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
777:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
778:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
779:                                                       void (*)(PetscInt, PetscInt, PetscInt,
780:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
781:                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
782:                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
783: PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
784:                                                     void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
785: PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
786:                                                     void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
787: PETSC_EXTERN PetscErrorCode PetscDSGetUpdate(PetscDS, PetscInt,
788:                                              void (**)(PetscInt, PetscInt, PetscInt,
789:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
790:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
791:                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
792: PETSC_EXTERN PetscErrorCode PetscDSSetUpdate(PetscDS, PetscInt,
793:                                              void (*)(PetscInt, PetscInt, PetscInt,
794:                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
795:                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
796:                                                       PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
797: PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void *);
798: PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
799: PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
800:                                                  void (**)(PetscInt, PetscInt, PetscInt,
801:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
802:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
803:                                                            PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
804:                                                  void (**)(PetscInt, PetscInt, PetscInt,
805:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
806:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
807:                                                            PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
808: PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
809:                                                  void (*)(PetscInt, PetscInt, PetscInt,
810:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
811:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
812:                                                           PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
813:                                                  void (*)(PetscInt, PetscInt, PetscInt,
814:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
815:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
816:                                                           PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
817: PETSC_EXTERN PetscErrorCode PetscDSHasBdJacobian(PetscDS, PetscBool *);
818: PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
819:                                                  void (**)(PetscInt, PetscInt, PetscInt,
820:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
821:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
822:                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
823:                                                  void (**)(PetscInt, PetscInt, PetscInt,
824:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
825:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
826:                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
827:                                                  void (**)(PetscInt, PetscInt, PetscInt,
828:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
829:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
830:                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
831:                                                  void (**)(PetscInt, PetscInt, PetscInt,
832:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
833:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
834:                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
835: PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
836:                                                  void (*)(PetscInt, PetscInt, PetscInt,
837:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
838:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
839:                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
840:                                                  void (*)(PetscInt, PetscInt, PetscInt,
841:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
842:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
843:                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
844:                                                  void (*)(PetscInt, PetscInt, PetscInt,
845:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
846:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
847:                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
848:                                                  void (*)(PetscInt, PetscInt, PetscInt,
849:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
850:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
851:                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
852: PETSC_EXTERN PetscErrorCode PetscDSHasBdJacobianPreconditioner(PetscDS, PetscBool *);
853: PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
854:                                                                void (**)(PetscInt, PetscInt, PetscInt,
855:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
856:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
857:                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
858:                                                                void (**)(PetscInt, PetscInt, PetscInt,
859:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
860:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
861:                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
862:                                                                void (**)(PetscInt, PetscInt, PetscInt,
863:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
864:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
865:                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
866:                                                                void (**)(PetscInt, PetscInt, PetscInt,
867:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
868:                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
869:                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
870: PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
871:                                                                void (*)(PetscInt, PetscInt, PetscInt,
872:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
873:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
874:                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
875:                                                                void (*)(PetscInt, PetscInt, PetscInt,
876:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
877:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
878:                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
879:                                                                void (*)(PetscInt, PetscInt, PetscInt,
880:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
881:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
882:                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
883:                                                                void (*)(PetscInt, PetscInt, PetscInt,
884:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
885:                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
886:                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
887: PETSC_EXTERN PetscErrorCode PetscDSGetExactSolution(PetscDS, PetscInt, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void **);
888: PETSC_EXTERN PetscErrorCode PetscDSSetExactSolution(PetscDS, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void *);
889: PETSC_EXTERN PetscErrorCode PetscDSGetExactSolutionTimeDerivative(PetscDS, PetscInt, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void **);
890: PETSC_EXTERN PetscErrorCode PetscDSSetExactSolutionTimeDerivative(PetscDS, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void *);
891: PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscTabulation *[]);
892: PETSC_EXTERN PetscErrorCode PetscDSGetFaceTabulation(PetscDS, PetscTabulation *[]);
893: PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
894: PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
895: PETSC_EXTERN PetscErrorCode PetscDSGetWorkspace(PetscDS, PetscReal **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
896: PETSC_EXTERN PetscErrorCode PetscDSCopyConstants(PetscDS, PetscDS);
897: PETSC_EXTERN PetscErrorCode PetscDSCopyExactSolutions(PetscDS, PetscDS);
898: PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS);
899: PETSC_EXTERN PetscErrorCode PetscDSSelectDiscretizations(PetscDS, PetscInt, const PetscInt[], PetscDS);
900: PETSC_EXTERN PetscErrorCode PetscDSSelectEquations(PetscDS, PetscInt, const PetscInt[], PetscDS);
901: PETSC_EXTERN PetscErrorCode PetscDSAddBoundary(PetscDS, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
902: PETSC_EXTERN PetscErrorCode PetscDSAddBoundaryByName(PetscDS, DMBoundaryConditionType, const char[], const char[], PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
903: PETSC_EXTERN PetscErrorCode PetscDSUpdateBoundary(PetscDS, PetscInt, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *);
904: PETSC_EXTERN PetscErrorCode PetscDSGetNumBoundary(PetscDS, PetscInt *);
905: PETSC_EXTERN PetscErrorCode PetscDSGetBoundary(PetscDS, PetscInt, PetscWeakForm *, DMBoundaryConditionType *, const char *[], DMLabel *, PetscInt *, const PetscInt *[], PetscInt *, PetscInt *, const PetscInt *[], void (**)(void), void (**)(void), void **);
906: PETSC_EXTERN PetscErrorCode PetscDSCopyBoundary(PetscDS, PetscInt, const PetscInt[], PetscDS);
907: PETSC_EXTERN PetscErrorCode PetscDSDestroyBoundary(PetscDS);

909: #endif