xref: /petsc/include/petscts.h (revision 09ee84387ed3895c141f9576cd10ed3bddfd46d1)
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 PETSC_EXTERN_CXX_BEGIN
9 
10 /*S
11      TS - Abstract PETSc object that manages all time-steppers (ODE integrators)
12 
13    Level: beginner
14 
15   Concepts: ODE solvers
16 
17 .seealso:  TSCreate(), TSSetType(), TSType, SNES, KSP, PC
18 S*/
19 typedef struct _p_TS* TS;
20 
21 /*E
22     TSType - String with the name of a PETSc TS method or the creation function
23        with an optional dynamic library name, for example
24        http://www.mcs.anl.gov/petsc/lib.a:mytscreate()
25 
26    Level: beginner
27 
28 .seealso: TSSetType(), TS
29 E*/
30 #define TSType char*
31 #define TSEULER           "euler"
32 #define TSBEULER          "beuler"
33 #define TSPSEUDO          "pseudo"
34 #define TSCN              "cn"
35 #define TSSUNDIALS        "sundials"
36 #define TSRK              "rk"
37 #define TSPYTHON          "python"
38 #define TSTHETA           "theta"
39 #define TSALPHA           "alpha"
40 #define TSGL              "gl"
41 #define TSSSP             "ssp"
42 #define TSARKIMEX         "arkimex"
43 
44 /*E
45     TSProblemType - Determines the type of problem this TS object is to be used to solve
46 
47    Level: beginner
48 
49 .seealso: TSCreate()
50 E*/
51 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType;
52 
53 typedef enum {
54   TS_CONVERGED_ITERATING      = 0,
55   TS_CONVERGED_TIME           = 1,
56   TS_CONVERGED_ITS            = 2,
57   TS_DIVERGED_NONLINEAR_SOLVE = -1,
58   TS_DIVERGED_STEP_REJECTED   = -2
59 } TSConvergedReason;
60 extern const char *const*TSConvergedReasons;
61 
62 /* Logging support */
63 extern PetscClassId  TS_CLASSID;
64 
65 extern PetscErrorCode   TSInitializePackage(const char[]);
66 
67 extern PetscErrorCode   TSCreate(MPI_Comm,TS*);
68 extern PetscErrorCode   TSDestroy(TS*);
69 
70 extern PetscErrorCode   TSSetProblemType(TS,TSProblemType);
71 extern PetscErrorCode   TSGetProblemType(TS,TSProblemType*);
72 extern PetscErrorCode   TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**));
73 extern PetscErrorCode   TSMonitorCancel(TS);
74 
75 extern PetscErrorCode   TSSetOptionsPrefix(TS,const char[]);
76 extern PetscErrorCode   TSAppendOptionsPrefix(TS,const char[]);
77 extern PetscErrorCode   TSGetOptionsPrefix(TS,const char *[]);
78 extern PetscErrorCode   TSSetFromOptions(TS);
79 extern PetscErrorCode   TSSetUp(TS);
80 extern PetscErrorCode   TSReset(TS);
81 
82 extern PetscErrorCode   TSSetSolution(TS,Vec);
83 extern PetscErrorCode   TSGetSolution(TS,Vec*);
84 
85 extern PetscErrorCode   TSSetDuration(TS,PetscInt,PetscReal);
86 extern PetscErrorCode   TSGetDuration(TS,PetscInt*,PetscReal*);
87 extern PetscErrorCode   TSSetExactFinalTime(TS,PetscBool);
88 
89 extern PetscErrorCode   TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*);
90 extern PetscErrorCode   TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*);
91 extern PetscErrorCode   TSMonitorSolutionCreate(TS,PetscViewer,void**);
92 extern PetscErrorCode   TSMonitorSolutionDestroy(void**);
93 extern PetscErrorCode   TSStep(TS);
94 extern PetscErrorCode   TSSolve(TS,Vec,PetscReal*);
95 extern PetscErrorCode   TSGetConvergedReason(TS,TSConvergedReason*);
96 
97 extern PetscErrorCode   TSSetInitialTimeStep(TS,PetscReal,PetscReal);
98 extern PetscErrorCode   TSGetTimeStep(TS,PetscReal*);
99 extern PetscErrorCode   TSGetTime(TS,PetscReal*);
100 extern PetscErrorCode   TSSetTime(TS,PetscReal);
101 extern PetscErrorCode   TSGetTimeStepNumber(TS,PetscInt*);
102 extern PetscErrorCode   TSSetTimeStep(TS,PetscReal);
103 
104 typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*);
105 typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
106 extern PetscErrorCode   TSSetRHSFunction(TS,Vec,TSRHSFunction,void*);
107 extern PetscErrorCode   TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**);
108 extern PetscErrorCode   TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*);
109 extern PetscErrorCode   TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**);
110 
111 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*);
112 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
113 extern PetscErrorCode   TSSetIFunction(TS,Vec,TSIFunction,void*);
114 extern PetscErrorCode   TSGetIFunction(TS,Vec*,TSIFunction*,void**);
115 extern PetscErrorCode   TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*);
116 extern PetscErrorCode   TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**);
117 
118 extern PetscErrorCode   TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*);
119 extern PetscErrorCode   TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
120 extern PetscErrorCode   TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*);
121 extern PetscErrorCode   TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
122 extern PetscErrorCode   TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
123 extern PetscErrorCode   TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
124 
125 extern PetscErrorCode   TSSetPreStep(TS, PetscErrorCode (*)(TS));
126 extern PetscErrorCode   TSSetPostStep(TS, PetscErrorCode (*)(TS));
127 extern PetscErrorCode   TSPreStep(TS);
128 extern PetscErrorCode   TSPostStep(TS);
129 extern PetscErrorCode   TSSetRetainStages(TS,PetscBool);
130 extern PetscErrorCode   TSInterpolate(TS,PetscReal,Vec);
131 
132 extern PetscErrorCode   TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*);
133 extern PetscErrorCode   TSPseudoDefaultTimeStep(TS,PetscReal*,void*);
134 extern PetscErrorCode   TSPseudoComputeTimeStep(TS,PetscReal *);
135 
136 extern PetscErrorCode   TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*);
137 extern PetscErrorCode   TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *);
138 extern PetscErrorCode   TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *);
139 extern PetscErrorCode   TSPseudoSetTimeStepIncrement(TS,PetscReal);
140 extern PetscErrorCode   TSPseudoIncrementDtFromInitialDt(TS);
141 
142 extern PetscErrorCode   TSPythonSetType(TS,const char[]);
143 
144 extern PetscErrorCode   TSComputeRHSFunction(TS,PetscReal,Vec,Vec);
145 extern PetscErrorCode   TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*);
146 extern PetscErrorCode   TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool);
147 extern PetscErrorCode   TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool);
148 
149 extern PetscErrorCode   TSVISetVariableBounds(TS,Vec,Vec);
150 
151 /* Dynamic creation and loading functions */
152 extern PetscFList TSList;
153 extern PetscBool  TSRegisterAllCalled;
154 extern PetscErrorCode   TSGetType(TS,const TSType*);
155 extern PetscErrorCode   TSSetType(TS,const TSType);
156 extern PetscErrorCode   TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS));
157 extern PetscErrorCode   TSRegisterAll(const char[]);
158 extern PetscErrorCode   TSRegisterDestroy(void);
159 
160 /*MC
161   TSRegisterDynamic - Adds a creation method to the TS package.
162 
163   Synopsis:
164   PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS))
165 
166   Not Collective
167 
168   Input Parameters:
169 + name        - The name of a new user-defined creation routine
170 . path        - The path (either absolute or relative) of the library containing this routine
171 . func_name   - The name of the creation routine
172 - create_func - The creation routine itself
173 
174   Notes:
175   TSRegisterDynamic() may be called multiple times to add several user-defined tses.
176 
177   If dynamic libraries are used, then the fourth input argument (create_func) is ignored.
178 
179   Sample usage:
180 .vb
181   TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate);
182 .ve
183 
184   Then, your ts type can be chosen with the procedural interface via
185 .vb
186     TS ts;
187     TSCreate(MPI_Comm, &ts);
188     TSSetType(ts, "my_ts")
189 .ve
190   or at runtime via the option
191 .vb
192     -ts_type my_ts
193 .ve
194 
195   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
196         If your function is not being put into a shared library then use TSRegister() instead
197 
198   Level: advanced
199 
200 .keywords: TS, register
201 .seealso: TSRegisterAll(), TSRegisterDestroy()
202 M*/
203 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
204 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0)
205 #else
206 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d)
207 #endif
208 
209 extern PetscErrorCode   TSGetSNES(TS,SNES*);
210 extern PetscErrorCode   TSGetKSP(TS,KSP*);
211 
212 extern PetscErrorCode   TSView(TS,PetscViewer);
213 
214 extern PetscErrorCode   TSSetApplicationContext(TS,void *);
215 extern PetscErrorCode   TSGetApplicationContext(TS,void *);
216 
217 extern PetscErrorCode   TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *);
218 extern PetscErrorCode   TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *);
219 extern PetscErrorCode   TSMonitorLGDestroy(PetscDrawLG*);
220 
221 /*E
222    TSSSPType - string with the name of TSSSP scheme.
223 
224    Level: beginner
225 
226 .seealso: TSSSPSetType(), TS
227 E*/
228 #define TSSSPType char*
229 #define TSSSPRKS2  "rks2"
230 #define TSSSPRKS3  "rks3"
231 #define TSSSPRK104 "rk104"
232 
233 extern PetscErrorCode TSSSPSetType(TS,const TSSSPType);
234 extern PetscErrorCode TSSSPGetType(TS,const TSSSPType*);
235 extern PetscErrorCode TSSSPSetNumStages(TS,PetscInt);
236 extern PetscErrorCode TSSSPGetNumStages(TS,PetscInt*);
237 
238 /*S
239    TSGLAdapt - Abstract object that manages time-step adaptivity
240 
241    Level: beginner
242 
243 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
244 S*/
245 typedef struct _p_TSGLAdapt *TSGLAdapt;
246 
247 /*E
248     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
249        with an optional dynamic library name, for example
250        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
251 
252    Level: beginner
253 
254 .seealso: TSGLAdaptSetType(), TS
255 E*/
256 #define TSGLAdaptType  char*
257 #define TSGLADAPT_NONE "none"
258 #define TSGLADAPT_SIZE "size"
259 #define TSGLADAPT_BOTH "both"
260 
261 /*MC
262    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
263 
264    Synopsis:
265    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
266 
267    Not Collective
268 
269    Input Parameters:
270 +  name_scheme - name of user-defined adaptivity scheme
271 .  path - path (either absolute or relative) the library containing this scheme
272 .  name_create - name of routine to create method context
273 -  routine_create - routine to create method context
274 
275    Notes:
276    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
277 
278    If dynamic libraries are used, then the fourth input argument (routine_create)
279    is ignored.
280 
281    Sample usage:
282 .vb
283    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
284                             "MySchemeCreate",MySchemeCreate);
285 .ve
286 
287    Then, your scheme can be chosen with the procedural interface via
288 $     TSGLAdaptSetType(ts,"my_scheme")
289    or at runtime via the option
290 $     -ts_adapt_type my_scheme
291 
292    Level: advanced
293 
294    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
295           and others of the form ${any_environmental_variable} occuring in pathname will be
296           replaced with appropriate values.
297 
298 .keywords: TSGLAdapt, register
299 
300 .seealso: TSGLAdaptRegisterAll()
301 M*/
302 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
303 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
304 #else
305 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
306 #endif
307 
308 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
309 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
310 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
311 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
312 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
313 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
314 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
315 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
316 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
317 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
318 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
319 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
320 
321 /*E
322     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
323        with an optional dynamic library name, for example
324        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
325 
326    Level: beginner
327 
328 .seealso: TSGLSetAcceptType(), TS
329 E*/
330 #define TSGLAcceptType  char*
331 #define TSGLACCEPT_ALWAYS "always"
332 
333 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
334 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
335 
336 /*MC
337    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
338 
339    Synopsis:
340    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
341 
342    Not Collective
343 
344    Input Parameters:
345 +  name_scheme - name of user-defined acceptance scheme
346 .  path - path (either absolute or relative) the library containing this scheme
347 .  name_create - name of routine to create method context
348 -  routine_create - routine to create method context
349 
350    Notes:
351    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
352 
353    If dynamic libraries are used, then the fourth input argument (routine_create)
354    is ignored.
355 
356    Sample usage:
357 .vb
358    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
359                              "MySchemeCreate",MySchemeCreate);
360 .ve
361 
362    Then, your scheme can be chosen with the procedural interface via
363 $     TSGLSetAcceptType(ts,"my_scheme")
364    or at runtime via the option
365 $     -ts_gl_accept_type my_scheme
366 
367    Level: advanced
368 
369    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
370           and others of the form ${any_environmental_variable} occuring in pathname will be
371           replaced with appropriate values.
372 
373 .keywords: TSGL, TSGLAcceptType, register
374 
375 .seealso: TSGLRegisterAll()
376 M*/
377 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
378 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
379 #else
380 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
381 #endif
382 
383 /*E
384   TSGLType - family of time integration method within the General Linear class
385 
386   Level: beginner
387 
388 .seealso: TSGLSetType(), TSGLRegister()
389 E*/
390 #define TSGLType char*
391 #define TSGL_IRKS   "irks"
392 
393 /*MC
394    TSGLRegisterDynamic - adds a TSGL implementation
395 
396    Synopsis:
397    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
398 
399    Not Collective
400 
401    Input Parameters:
402 +  name_scheme - name of user-defined general linear scheme
403 .  path - path (either absolute or relative) the library containing this scheme
404 .  name_create - name of routine to create method context
405 -  routine_create - routine to create method context
406 
407    Notes:
408    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
409 
410    If dynamic libraries are used, then the fourth input argument (routine_create)
411    is ignored.
412 
413    Sample usage:
414 .vb
415    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
416                        "MySchemeCreate",MySchemeCreate);
417 .ve
418 
419    Then, your scheme can be chosen with the procedural interface via
420 $     TSGLSetType(ts,"my_scheme")
421    or at runtime via the option
422 $     -ts_gl_type my_scheme
423 
424    Level: advanced
425 
426    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
427           and others of the form ${any_environmental_variable} occuring in pathname will be
428           replaced with appropriate values.
429 
430 .keywords: TSGL, register
431 
432 .seealso: TSGLRegisterAll()
433 M*/
434 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
435 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
436 #else
437 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
438 #endif
439 
440 extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
441 extern PetscErrorCode  TSGLRegisterAll(const char[]);
442 extern PetscErrorCode  TSGLRegisterDestroy(void);
443 extern PetscErrorCode  TSGLInitializePackage(const char[]);
444 extern PetscErrorCode  TSGLFinalizePackage(void);
445 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
446 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
447 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
448 
449 #define TSARKIMEXType char*
450 #define TSARKIMEX2D "2d"
451 #define TSARKIMEX2E "2e"
452 #define TSARKIMEX3  "3"
453 #define TSARKIMEX4  "4"
454 #define TSARKIMEX5  "5"
455 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
456 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
457 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
458 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
459 extern PetscErrorCode TSARKIMEXFinalizePackage(void);
460 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
461 extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
462 extern PetscErrorCode TSARKIMEXRegisterAll(void);
463 
464 /*
465        PETSc interface to Sundials
466 */
467 #ifdef PETSC_HAVE_SUNDIALS
468 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
469 extern const char *TSSundialsLmmTypes[];
470 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
471 extern const char *TSSundialsGramSchmidtTypes[];
472 extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
473 extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
474 extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
475 extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
476 extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
477 extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
478 extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
479 extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
480 extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
481 extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
482 extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
483 extern PetscErrorCode   TSSundialsSetMaxl(TS,PetscInt);
484 #endif
485 
486 extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
487 
488 extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
489 extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
490 extern PetscErrorCode  TSThetaGetEndpoint(TS,PetscBool*);
491 extern PetscErrorCode  TSThetaSetEndpoint(TS,PetscBool);
492 
493 extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
494 extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
495 extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
496 extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
497 extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
498 
499 extern PetscErrorCode  TSSetDM(TS,DM);
500 extern PetscErrorCode  TSGetDM(TS,DM*);
501 
502 extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
503 extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
504 
505 PETSC_EXTERN_CXX_END
506 #endif
507