xref: /petsc/include/petscts.h (revision 8efa2a282cafb3126d64e7372c7c93ad3b71bda1)
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 /*S
222    TSGLAdapt - Abstract object that manages time-step adaptivity
223 
224    Level: beginner
225 
226 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType
227 S*/
228 typedef struct _p_TSGLAdapt *TSGLAdapt;
229 
230 /*E
231     TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function
232        with an optional dynamic library name, for example
233        http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate()
234 
235    Level: beginner
236 
237 .seealso: TSGLAdaptSetType(), TS
238 E*/
239 #define TSGLAdaptType  char*
240 #define TSGLADAPT_NONE "none"
241 #define TSGLADAPT_SIZE "size"
242 #define TSGLADAPT_BOTH "both"
243 
244 /*MC
245    TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation
246 
247    Synopsis:
248    PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
249 
250    Not Collective
251 
252    Input Parameters:
253 +  name_scheme - name of user-defined adaptivity scheme
254 .  path - path (either absolute or relative) the library containing this scheme
255 .  name_create - name of routine to create method context
256 -  routine_create - routine to create method context
257 
258    Notes:
259    TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families.
260 
261    If dynamic libraries are used, then the fourth input argument (routine_create)
262    is ignored.
263 
264    Sample usage:
265 .vb
266    TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
267                             "MySchemeCreate",MySchemeCreate);
268 .ve
269 
270    Then, your scheme can be chosen with the procedural interface via
271 $     TSGLAdaptSetType(ts,"my_scheme")
272    or at runtime via the option
273 $     -ts_adapt_type my_scheme
274 
275    Level: advanced
276 
277    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
278           and others of the form ${any_environmental_variable} occuring in pathname will be
279           replaced with appropriate values.
280 
281 .keywords: TSGLAdapt, register
282 
283 .seealso: TSGLAdaptRegisterAll()
284 M*/
285 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
286 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,0)
287 #else
288 #  define TSGLAdaptRegisterDynamic(a,b,c,d)  TSGLAdaptRegister(a,b,c,d)
289 #endif
290 
291 extern PetscErrorCode  TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt));
292 extern PetscErrorCode  TSGLAdaptRegisterAll(const char[]);
293 extern PetscErrorCode  TSGLAdaptRegisterDestroy(void);
294 extern PetscErrorCode  TSGLAdaptInitializePackage(const char[]);
295 extern PetscErrorCode  TSGLAdaptFinalizePackage(void);
296 extern PetscErrorCode  TSGLAdaptCreate(MPI_Comm,TSGLAdapt*);
297 extern PetscErrorCode  TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType);
298 extern PetscErrorCode  TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]);
299 extern PetscErrorCode  TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *);
300 extern PetscErrorCode  TSGLAdaptView(TSGLAdapt,PetscViewer);
301 extern PetscErrorCode  TSGLAdaptSetFromOptions(TSGLAdapt);
302 extern PetscErrorCode  TSGLAdaptDestroy(TSGLAdapt*);
303 
304 /*E
305     TSGLAcceptType - String with the name of TSGLAccept scheme or the function
306        with an optional dynamic library name, for example
307        http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept()
308 
309    Level: beginner
310 
311 .seealso: TSGLSetAcceptType(), TS
312 E*/
313 #define TSGLAcceptType  char*
314 #define TSGLACCEPT_ALWAYS "always"
315 
316 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *);
317 extern PetscErrorCode  TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction);
318 
319 /*MC
320    TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme
321 
322    Synopsis:
323    PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
324 
325    Not Collective
326 
327    Input Parameters:
328 +  name_scheme - name of user-defined acceptance scheme
329 .  path - path (either absolute or relative) the library containing this scheme
330 .  name_create - name of routine to create method context
331 -  routine_create - routine to create method context
332 
333    Notes:
334    TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families.
335 
336    If dynamic libraries are used, then the fourth input argument (routine_create)
337    is ignored.
338 
339    Sample usage:
340 .vb
341    TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
342                              "MySchemeCreate",MySchemeCreate);
343 .ve
344 
345    Then, your scheme can be chosen with the procedural interface via
346 $     TSGLSetAcceptType(ts,"my_scheme")
347    or at runtime via the option
348 $     -ts_gl_accept_type my_scheme
349 
350    Level: advanced
351 
352    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
353           and others of the form ${any_environmental_variable} occuring in pathname will be
354           replaced with appropriate values.
355 
356 .keywords: TSGL, TSGLAcceptType, register
357 
358 .seealso: TSGLRegisterAll()
359 M*/
360 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
361 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0)
362 #else
363 #  define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d)
364 #endif
365 
366 /*E
367   TSGLType - family of time integration method within the General Linear class
368 
369   Level: beginner
370 
371 .seealso: TSGLSetType(), TSGLRegister()
372 E*/
373 #define TSGLType char*
374 #define TSGL_IRKS   "irks"
375 
376 /*MC
377    TSGLRegisterDynamic - adds a TSGL implementation
378 
379    Synopsis:
380    PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS))
381 
382    Not Collective
383 
384    Input Parameters:
385 +  name_scheme - name of user-defined general linear scheme
386 .  path - path (either absolute or relative) the library containing this scheme
387 .  name_create - name of routine to create method context
388 -  routine_create - routine to create method context
389 
390    Notes:
391    TSGLRegisterDynamic() may be called multiple times to add several user-defined families.
392 
393    If dynamic libraries are used, then the fourth input argument (routine_create)
394    is ignored.
395 
396    Sample usage:
397 .vb
398    TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a,
399                        "MySchemeCreate",MySchemeCreate);
400 .ve
401 
402    Then, your scheme can be chosen with the procedural interface via
403 $     TSGLSetType(ts,"my_scheme")
404    or at runtime via the option
405 $     -ts_gl_type my_scheme
406 
407    Level: advanced
408 
409    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
410           and others of the form ${any_environmental_variable} occuring in pathname will be
411           replaced with appropriate values.
412 
413 .keywords: TSGL, register
414 
415 .seealso: TSGLRegisterAll()
416 M*/
417 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
418 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,0)
419 #else
420 #  define TSGLRegisterDynamic(a,b,c,d)       TSGLRegister(a,b,c,d)
421 #endif
422 
423 extern PetscErrorCode  TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS));
424 extern PetscErrorCode  TSGLRegisterAll(const char[]);
425 extern PetscErrorCode  TSGLRegisterDestroy(void);
426 extern PetscErrorCode  TSGLInitializePackage(const char[]);
427 extern PetscErrorCode  TSGLFinalizePackage(void);
428 extern PetscErrorCode  TSGLSetType(TS,const TSGLType);
429 extern PetscErrorCode  TSGLGetAdapt(TS,TSGLAdapt*);
430 extern PetscErrorCode  TSGLSetAcceptType(TS,const TSGLAcceptType);
431 
432 #define TSARKIMEXType char*
433 #define TSARKIMEX2D "2d"
434 #define TSARKIMEX2E "2e"
435 #define TSARKIMEX3  "3"
436 #define TSARKIMEX4  "4"
437 #define TSARKIMEX5  "5"
438 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*);
439 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType);
440 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool);
441 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]);
442 extern PetscErrorCode TSARKIMEXFinalizePackage(void);
443 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]);
444 extern PetscErrorCode TSARKIMEXRegisterDestroy(void);
445 extern PetscErrorCode TSARKIMEXRegisterAll(void);
446 
447 /*
448        PETSc interface to Sundials
449 */
450 #ifdef PETSC_HAVE_SUNDIALS
451 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType;
452 extern const char *TSSundialsLmmTypes[];
453 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType;
454 extern const char *TSSundialsGramSchmidtTypes[];
455 extern PetscErrorCode   TSSundialsSetType(TS,TSSundialsLmmType);
456 extern PetscErrorCode   TSSundialsGetPC(TS,PC*);
457 extern PetscErrorCode   TSSundialsSetTolerance(TS,PetscReal,PetscReal);
458 extern PetscErrorCode   TSSundialsSetMinTimeStep(TS,PetscReal);
459 extern PetscErrorCode   TSSundialsSetMaxTimeStep(TS,PetscReal);
460 extern PetscErrorCode   TSSundialsGetIterations(TS,PetscInt *,PetscInt *);
461 extern PetscErrorCode   TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType);
462 extern PetscErrorCode   TSSundialsSetGMRESRestart(TS,PetscInt);
463 extern PetscErrorCode   TSSundialsSetLinearTolerance(TS,PetscReal);
464 extern PetscErrorCode   TSSundialsMonitorInternalSteps(TS,PetscBool );
465 extern PetscErrorCode   TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]);
466 extern PetscErrorCode   TSSundialsSetMaxl(TS,PetscInt);
467 #endif
468 
469 extern PetscErrorCode   TSRKSetTolerance(TS,PetscReal);
470 
471 extern PetscErrorCode  TSThetaSetTheta(TS,PetscReal);
472 extern PetscErrorCode  TSThetaGetTheta(TS,PetscReal*);
473 extern PetscErrorCode  TSThetaGetEndpoint(TS,PetscBool*);
474 extern PetscErrorCode  TSThetaSetEndpoint(TS,PetscBool);
475 
476 extern PetscErrorCode  TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*);
477 extern PetscErrorCode  TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*);
478 extern PetscErrorCode  TSAlphaSetRadius(TS,PetscReal);
479 extern PetscErrorCode  TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal);
480 extern PetscErrorCode  TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*);
481 
482 extern PetscErrorCode  TSSetDM(TS,DM);
483 extern PetscErrorCode  TSGetDM(TS,DM*);
484 
485 extern PetscErrorCode  SNESTSFormFunction(SNES,Vec,Vec,void*);
486 extern PetscErrorCode  SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
487 
488 PETSC_EXTERN_CXX_END
489 #endif
490