xref: /petsc/include/petscts.h (revision abadfa0ff52f084808973d56e72ca5b14b37a6b3)
1 /*
2    User interface for the timestepping package. This package
3    is for use in solving time-dependent PDEs.
4 */
5 #if !defined(__PETSCTS_H)
6 #define __PETSCTS_H
7 #include <petscsnes.h>
8 
9 /*S
10      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
11 
12    Level: beginner
13 
14   Concepts: ODE solvers
15 
16 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
17 S*/
18 typedef struct _p_TS* TS;
19 
20 /*J
21     TSType - String with the name of a PETSc TS method.
22 
23    Level: beginner
24 
25 .seealso: TSSetType(), TS, TSRegister()
26 J*/
27 typedef const char* TSType;
28 #define TSEULER           "euler"
29 #define TSBEULER          "beuler"
30 #define TSPSEUDO          "pseudo"
31 #define TSCN              "cn"
32 #define TSSUNDIALS        "sundials"
33 #define TSRK              "rk"
34 #define TSPYTHON          "python"
35 #define TSTHETA           "theta"
36 #define TSALPHA           "alpha"
37 #define TSGL              "gl"
38 #define TSSSP             "ssp"
39 #define TSARKIMEX         "arkimex"
40 #define TSROSW            "rosw"
41 #define TSEIMEX           "eimex"
42 #define TSMIMEX           "mimex"
43 /*E
44     TSProblemType - Determines the type of problem this TS object is to be used to solve
45 
46    Level: beginner
47 
48 .seealso: TSCreate()
49 E*/
50 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
51 
52 /*E
53    TSEquationType - type of TS problem that is solved
54 
55    Level: beginner
56 
57    Developer Notes: this must match petsc-finclude/petscts.h
58 
59    Supported types are:
60      TS_EQ_UNSPECIFIED (default)
61      TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0
62      TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0
63 
64 .seealso: TSGetEquationType(), TSSetEquationType()
65 E*/
66 typedef enum {
67   TS_EQ_UNSPECIFIED               = -1,
68   TS_EQ_EXPLICIT                  = 0,
69   TS_EQ_ODE_EXPLICIT              = 1,
70   TS_EQ_DAE_SEMI_EXPLICIT_INDEX1  = 100,
71   TS_EQ_DAE_SEMI_EXPLICIT_INDEX2  = 200,
72   TS_EQ_DAE_SEMI_EXPLICIT_INDEX3  = 300,
73   TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500,
74   TS_EQ_IMPLICIT                  = 1000,
75   TS_EQ_ODE_IMPLICIT              = 1001,
76   TS_EQ_DAE_IMPLICIT_INDEX1       = 1100,
77   TS_EQ_DAE_IMPLICIT_INDEX2       = 1200,
78   TS_EQ_DAE_IMPLICIT_INDEX3       = 1300,
79   TS_EQ_DAE_IMPLICIT_INDEXHI      = 1500
80 } TSEquationType;
81 PETSC_EXTERN const char *const*TSEquationTypes;
82 
83 /*E
84    TSConvergedReason - reason a TS method has converged or not
85 
86    Level: beginner
87 
88    Developer Notes: this must match petsc-finclude/petscts.h
89 
90    Each reason has its own manual page.
91 
92 .seealso: TSGetConvergedReason()
93 E*/
94 typedef enum {
95   TS_CONVERGED_ITERATING      = 0,
96   TS_CONVERGED_TIME           = 1,
97   TS_CONVERGED_ITS            = 2,
98   TS_CONVERGED_USER           = 3,
99   TS_CONVERGED_EVENT          = 4,
100   TS_DIVERGED_NONLINEAR_SOLVE = -1,
101   TS_DIVERGED_STEP_REJECTED   = -2
102 } TSConvergedReason;
103 PETSC_EXTERN const char *const*TSConvergedReasons;
104 
105 /*MC
106    TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve()
107 
108    Level: beginner
109 
110 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt()
111 M*/
112 
113 /*MC
114    TS_CONVERGED_TIME - the final time was reached
115 
116    Level: beginner
117 
118 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime()
119 M*/
120 
121 /*MC
122    TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time
123 
124    Level: beginner
125 
126 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration()
127 M*/
128 
129 /*MC
130    TS_CONVERGED_USER - user requested termination
131 
132    Level: beginner
133 
134 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
135 M*/
136 
137 /*MC
138    TS_CONVERGED_EVENT - user requested termination on event detection
139 
140    Level: beginner
141 
142 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration()
143 M*/
144 
145 /*MC
146    TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed
147 
148    Level: beginner
149 
150    Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures.
151 
152 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures()
153 M*/
154 
155 /*MC
156    TS_DIVERGED_STEP_REJECTED - too many steps were rejected
157 
158    Level: beginner
159 
160    Notes: See TSSetMaxStepRejections() for how to allow more step rejections.
161 
162 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections()
163 M*/
164 
165 /*E
166    TSExactFinalTimeOption - option for handling of final time step
167 
168    Level: beginner
169 
170    Developer Notes: this must match petsc-finclude/petscts.h
171 
172 $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
173 $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
174 $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
175 .seealso: TSGetConvergedReason(), TSSetExactFinalTime()
176 
177 E*/
178 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption;
179 PETSC_EXTERN const char *const TSExactFinalTimeOptions[];
180 
181 
182 /* Logging support */
183 PETSC_EXTERN PetscClassId TS_CLASSID;
184 PETSC_EXTERN PetscClassId DMTS_CLASSID;
185 
186 PETSC_EXTERN PetscErrorCode TSInitializePackage(void);
187 
188 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*);
189 PETSC_EXTERN PetscErrorCode TSDestroy(TS*);
190 
191 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType);
192 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*);
193 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec);
194 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
195 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS);
196 
197 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]);
198 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]);
199 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]);
200 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS);
201 PETSC_EXTERN PetscErrorCode TSSetUp(TS);
202 PETSC_EXTERN PetscErrorCode TSReset(TS);
203 
204 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec);
205 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*);
206 
207 /*S
208      TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage)
209 
210    Level: advanced
211 
212   Concepts: ODE solvers, adjoint
213 
214 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy()
215 S*/
216 typedef struct _p_TSTrajectory* TSTrajectory;
217 
218 /*J
219     TSTrajectoryType - String with the name of a PETSc TS trajectory storage method
220 
221    Level: intermediate
222 
223 .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType()
224 J*/
225 typedef const char* TSTrajectoryType;
226 #define TSTRAJECTORYBASIC      "basic"
227 
228 PETSC_EXTERN PetscFunctionList TSTrajectoryList;
229 PETSC_EXTERN PetscClassId      TSTRAJECTORY_CLASSID;
230 PETSC_EXTERN PetscBool         TSTrajectoryRegisterAllCalled;
231 
232 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS);
233 
234 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*);
235 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*);
236 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,const TSTrajectoryType);
237 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec);
238 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal);
239 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory);
240 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void);
241 
242 PETSC_EXTERN PetscErrorCode TSAdjointSetGradients(TS,PetscInt,Vec*,Vec*);
243 PETSC_EXTERN PetscErrorCode TSAdjointGetGradients(TS,PetscInt*,Vec**,Vec**);
244 
245 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*);
246 PETSC_EXTERN PetscErrorCode TSAdjointSetCostIntegrand(TS,PetscInt,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),Vec*,PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),Vec*,PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),void*);
247 PETSC_EXTERN PetscErrorCode TSAdjointGetCostIntegral(TS,Vec*);
248 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS);
249 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt);
250 
251 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat);
252 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS);
253 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS);
254 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*);
255 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*);
256 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec);
257 
258 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal);
259 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*);
260 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption);
261 
262 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
263 
264 typedef struct _n_TSMonitorDrawCtx*  TSMonitorDrawCtx;
265 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *);
266 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*);
267 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*);
268 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*);
269 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*);
270 
271 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*);
272 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*);
273 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*);
274 
275 PETSC_EXTERN PetscErrorCode TSStep(TS);
276 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*);
277 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec);
278 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*);
279 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType);
280 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*);
281 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason);
282 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*);
283 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*);
284 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*);
285 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*);
286 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt);
287 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*);
288 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt);
289 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool);
290 PETSC_EXTERN PetscErrorCode TSRollBack(TS);
291 PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*);
292 
293 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**);
294 
295 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal);
296 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*);
297 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*);
298 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal);
299 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*);
300 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal);
301 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*);
302 
303 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
304 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*);
305 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
306 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
307 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
308 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
309 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool);
310 
311 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*);
312 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*);
313 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
314 
315 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
316 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
317 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*);
318 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**);
319 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
320 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
321 
322 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
323 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*);
324 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
325 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
326 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec);
327 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec);
328 
329 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS));
330 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal));
331 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*));
332 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS));
333 PETSC_EXTERN PetscErrorCode TSPreStep(TS);
334 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal);
335 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*);
336 PETSC_EXTERN PetscErrorCode TSPostStep(TS);
337 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool);
338 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec);
339 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec);
340 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*);
341 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*);
342 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal);
343 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*);
344 
345 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
346 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*);
347 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *);
348 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal);
349 
350 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
351 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *);
352 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
353 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal);
354 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS);
355 
356 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]);
357 
358 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
359 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat);
360 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
361 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool);
362 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
363 
364 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec);
365 
366 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*);
367 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**);
368 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*);
369 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**);
370 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*);
371 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**);
372 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*);
373 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**);
374 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*);
375 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**);
376 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*);
377 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**);
378 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*);
379 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal);
380 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, void (**)(const PetscReal[], PetscScalar *, void *), void **);
381 
382 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*);
383 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*);
384 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*);
385 
386 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
387 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer));
388 
389 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
390 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
391 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
392 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
393 
394 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *);
395 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *);
396 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *);
397 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *);
398 
399 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
400 
401 typedef struct _n_TSMonitorLGCtx*  TSMonitorLGCtx;
402 typedef struct {
403   Vec            ray;
404   VecScatter     scatter;
405   PetscViewer    viewer;
406   TSMonitorLGCtx lgctx;
407 } TSMonitorDMDARayCtx;
408 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**);
409 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*);
410 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*);
411 
412 
413 /* Dynamic creation and loading functions */
414 PETSC_EXTERN PetscFunctionList TSList;
415 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*);
416 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType);
417 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS));
418 
419 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*);
420 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES);
421 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*);
422 
423 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer);
424 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer);
425 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
426 
427 
428 #define TS_FILE_CLASSID 1211225
429 
430 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *);
431 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *);
432 
433 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *);
434 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*);
435 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *);
436 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *);
437 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *);
438 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *);
439 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *);
440 
441 typedef struct _n_TSMonitorSPEigCtx*  TSMonitorSPEigCtx;
442 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *);
443 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*);
444 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *);
445 
446 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,void*),void*);
447 /*J
448    TSSSPType - string with the name of TSSSP scheme.
449 
450    Level: beginner
451 
452 .seealso: TSSSPSetType(), TS
453 J*/
454 typedef const char* TSSSPType;
455 #define TSSSPRKS2  "rks2"
456 #define TSSSPRKS3  "rks3"
457 #define TSSSPRK104 "rk104"
458 
459 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType);
460 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*);
461 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
462 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
463 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void);
464 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void);
465 
466 /*S
467    TSAdapt - Abstract object that manages time-step adaptivity
468 
469    Level: beginner
470 
471 .seealso: TS, TSAdaptCreate(), TSAdaptType
472 S*/
473 typedef struct _p_TSAdapt *TSAdapt;
474 
475 /*E
476     TSAdaptType - String with the name of TSAdapt scheme.
477 
478    Level: beginner
479 
480 .seealso: TSAdaptSetType(), TS
481 E*/
482 typedef const char *TSAdaptType;
483 #define TSADAPTBASIC "basic"
484 #define TSADAPTNONE  "none"
485 #define TSADAPTCFL   "cfl"
486 
487 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*);
488 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt));
489 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void);
490 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void);
491 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*);
492 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType);
493 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]);
494 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt);
495 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool);
496 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**);
497 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*);
498 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*);
499 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer);
500 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer);
501 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt);
502 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*);
503 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool);
504 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal);
505 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*));
506 
507 /*S
508    TSGLAdapt - Abstract object that manages time-step adaptivity
509 
510    Level: beginner
511 
512    Developer Notes:
513    This functionality should be replaced by the TSAdapt.
514 
515 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
516 S*/
517 typedef struct _p_TSGLAdapt *TSGLAdapt;
518 
519 /*J
520     TSGLAdaptType - String with the name of TSGLAdapt scheme
521 
522    Level: beginner
523 
524 .seealso: TSGLAdaptSetType(), TS
525 J*/
526 typedef const char *TSGLAdaptType;
527 #define TSGLADAPT_NONE "none"
528 #define TSGLADAPT_SIZE "size"
529 #define TSGLADAPT_BOTH "both"
530 
531 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt));
532 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void);
533 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void);
534 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
535 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType);
536 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
537 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
538 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer);
539 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt);
540 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*);
541 
542 /*J
543     TSGLAcceptType - String with the name of TSGLAccept scheme
544 
545    Level: beginner
546 
547 .seealso: TSGLSetAcceptType(), TS
548 J*/
549 typedef const char *TSGLAcceptType;
550 #define TSGLACCEPT_ALWAYS "always"
551 
552 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
553 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction);
554 
555 /*J
556   TSGLType - family of time integration method within the General Linear class
557 
558   Level: beginner
559 
560 .seealso: TSGLSetType(), TSGLRegister()
561 J*/
562 typedef const char* TSGLType;
563 #define TSGL_IRKS   "irks"
564 
565 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS));
566 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void);
567 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void);
568 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType);
569 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*);
570 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType);
571 
572 /*J
573     TSEIMEXType - String with the name of an Extrapolated IMEX method.
574 
575    Level: beginner
576 
577 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister()
578 J*/
579 #define TSEIMEXType   char*
580 
581 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt);
582 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt);
583 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool);
584 
585 /*J
586     TSRKType - String with the name of a Runge-Kutta method.
587 
588    Level: beginner
589 
590 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister()
591 J*/
592 typedef const char* TSRKType;
593 #define TSRK1FE   "1fe"
594 #define TSRK2A    "2a"
595 #define TSRK3     "3"
596 #define TSRK3BS   "3bs"
597 #define TSRK4     "4"
598 #define TSRK5F    "5f"
599 #define TSRK5DP   "5dp"
600 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*);
601 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType);
602 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool);
603 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
604 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void);
605 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void);
606 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void);
607 
608 /*J
609     TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method.
610 
611    Level: beginner
612 
613 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister()
614 J*/
615 typedef const char* TSARKIMEXType;
616 #define TSARKIMEX1BEE   "1bee"
617 #define TSARKIMEXA2     "a2"
618 #define TSARKIMEXL2     "l2"
619 #define TSARKIMEXARS122 "ars122"
620 #define TSARKIMEX2C     "2c"
621 #define TSARKIMEX2D     "2d"
622 #define TSARKIMEX2E     "2e"
623 #define TSARKIMEXPRSSP2 "prssp2"
624 #define TSARKIMEX3      "3"
625 #define TSARKIMEXBPR3   "bpr3"
626 #define TSARKIMEXARS443 "ars443"
627 #define TSARKIMEX4      "4"
628 #define TSARKIMEX5      "5"
629 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*);
630 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType);
631 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
632 PETSC_EXTERN PetscErrorCode TSARKIMEXRegister(TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
633 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void);
634 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void);
635 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void);
636 
637 /*J
638     TSRosWType - String with the name of a Rosenbrock-W method.
639 
640    Level: beginner
641 
642 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister()
643 J*/
644 typedef const char* TSRosWType;
645 #define TSROSW2M          "2m"
646 #define TSROSW2P          "2p"
647 #define TSROSWRA3PW       "ra3pw"
648 #define TSROSWRA34PW2     "ra34pw2"
649 #define TSROSWRODAS3      "rodas3"
650 #define TSROSWSANDU3      "sandu3"
651 #define TSROSWASSP3P3S1C  "assp3p3s1c"
652 #define TSROSWLASSP3P4S2C "lassp3p4s2c"
653 #define TSROSWLLSSP3P4S2C "llssp3p4s2c"
654 #define TSROSWARK3        "ark3"
655 #define TSROSWTHETA1      "theta1"
656 #define TSROSWTHETA2      "theta2"
657 #define TSROSWGRK4T       "grk4t"
658 #define TSROSWSHAMP4      "shamp4"
659 #define TSROSWVELDD4      "veldd4"
660 #define TSROSW4L          "4l"
661 
662 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*);
663 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType);
664 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool);
665 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]);
666 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
667 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void);
668 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void);
669 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void);
670 
671 /*
672        PETSc interface to Sundials
673 */
674 #ifdef PETSC_HAVE_SUNDIALS
675 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
676 PETSC_EXTERN const char *const TSSundialsLmmTypes[];
677 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
678 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[];
679 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType);
680 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*);
681 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal);
682 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal);
683 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal);
684 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
685 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
686 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt);
687 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal);
688 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool );
689 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
690 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt);
691 #endif
692 
693 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal);
694 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*);
695 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*);
696 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool);
697 
698 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
699 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
700 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal);
701 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
702 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
703 
704 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM);
705 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*);
706 
707 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*);
708 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*);
709 
710 #endif
711