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